xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 8d298b4ff90fbe92523c3b0af19bac1dc5de81d2)
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*/
11e9e4536cSHong Zhang /*
12e9e4536cSHong Zhang #define DEBUG_MATMATMULT
13e9e4536cSHong Zhang  */
14ec01deb9SMatthew Knepley EXTERN_C_BEGIN
156284ec50SHong Zhang #undef __FUNCT__
165c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1738baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1838baddfdSBarry Smith {
19dfbe8321SBarry Smith   PetscErrorCode ierr;
2025023028SHong Zhang   PetscBool      scalable=PETSC_FALSE;
215c66b693SKris Buschelman 
225c66b693SKris Buschelman   PetscFunctionBegin;
2326be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
24598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2525023028SHong Zhang     ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_scalable",&scalable,PETSC_NULL);CHKERRQ(ierr);
2625023028SHong Zhang     if (scalable){
2725023028SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
2825023028SHong Zhang     } else {
2926be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3025023028SHong Zhang     }
31598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3226be0446SHong Zhang   }
33598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3401e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
35598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
365c66b693SKris Buschelman   PetscFunctionReturn(0);
375c66b693SKris Buschelman }
38ec01deb9SMatthew Knepley EXTERN_C_END
391c24bd37SHong Zhang 
40e9e4536cSHong Zhang /*
41e9e4536cSHong Zhang  MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ - Get symbolic structure of C=A*B
42e9e4536cSHong Zhang   Input Parameter:
43e9e4536cSHong Zhang .    am, Ai, Aj - number of rows and structure of A
44e9e4536cSHong Zhang .    bm, bn, Bi, Bj - number of rows, columns, and structure of B
45e9e4536cSHong Zhang .    fill - filll ratio See MatMatMult()
46e9e4536cSHong Zhang 
47e9e4536cSHong Zhang   Output Parameter:
48e9e4536cSHong Zhang .    Ci, Cj - structure of C = A*B
49e9e4536cSHong Zhang .    nspacedouble - number of extra mallocs
50e9e4536cSHong Zhang  */
5126be0446SHong Zhang #undef __FUNCT__
52b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ"
53971236abSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
5458c24d83SHong Zhang {
55dfbe8321SBarry Smith   PetscErrorCode     ierr;
56a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
57971236abSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data;
58971236abSHong Zhang   const PetscInt     *Ai=a->i,*Aj=a->j,*Bi=b->i,*Bj=b->j,am=A->rmap->N,bm=B->rmap->N,bn=B->cmap->N;
59971236abSHong Zhang   const PetscInt     *aj=Aj,*bi=Bi,*bj=Bj,*bjj;
60971236abSHong Zhang   PetscInt           *ci,*cj;
61b561aa9dSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0;
62be0fcf8dSHong Zhang   PetscBT            lnkbt;
6358c24d83SHong Zhang 
6458c24d83SHong Zhang   PetscFunctionBegin;
6558c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
6658c24d83SHong Zhang   /* free space for accumulating nonzero column info */
6738baddfdSBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6858c24d83SHong Zhang   ci[0] = 0;
6958c24d83SHong Zhang 
70be0fcf8dSHong Zhang   /* create and initialize a linked list */
71be0fcf8dSHong Zhang   nlnk = bn+1;
72be0fcf8dSHong Zhang   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
7358c24d83SHong Zhang 
74c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
75971236abSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(Ai[am]+Bi[bm])),&free_space);CHKERRQ(ierr);
7658c24d83SHong Zhang   current_space = free_space;
7758c24d83SHong Zhang 
78bf9555e6SHong Zhang   /* Determine ci and cj */
7958c24d83SHong Zhang   for (i=0; i<am; i++) {
80971236abSHong Zhang     anzi = Ai[i+1] - Ai[i];
8158c24d83SHong Zhang     cnzi = 0;
82971236abSHong Zhang     aj   = Aj + Ai[i];
8377584fe6SHong Zhang     for (j=0; j<anzi; j++){
84b561aa9dSHong Zhang       brow = aj[j];
8558c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
8658c24d83SHong Zhang       bjj  = bj + bi[brow];
871c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
88dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
891c239cc6SHong Zhang       cnzi += nlnk;
9058c24d83SHong Zhang     }
9158c24d83SHong Zhang 
9258c24d83SHong Zhang     /* If free space is not available, make more free space */
9358c24d83SHong Zhang     /* Double the amount of total space in the list */
9458c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
954238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
96b561aa9dSHong Zhang       ndouble++;
9758c24d83SHong Zhang     }
9858c24d83SHong Zhang 
99c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
100be0fcf8dSHong Zhang     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
101c5db241fSHong Zhang     current_space->array           += cnzi;
10258c24d83SHong Zhang     current_space->local_used      += cnzi;
10358c24d83SHong Zhang     current_space->local_remaining -= cnzi;
10458c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
10558c24d83SHong Zhang   }
10658c24d83SHong Zhang 
10758c24d83SHong Zhang   /* Column indices are in the list of free space */
10858c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
10958c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
11038baddfdSBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
111a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
112be0fcf8dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
11358c24d83SHong Zhang 
114b561aa9dSHong Zhang   *Ci           = ci;
115b561aa9dSHong Zhang   *Cj           = cj;
116b561aa9dSHong Zhang   *nspacedouble = ndouble;
117b561aa9dSHong Zhang   PetscFunctionReturn(0);
118b561aa9dSHong Zhang }
119b561aa9dSHong Zhang 
12031e66dd0SHong Zhang /*
12125023028SHong Zhang  MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable - same as MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ()
122bf9555e6SHong Zhang    but replaces O(bn) array 'lnkbt' with a scalable array 'abj' of size Crmax, or with a condensed list.
12331e66dd0SHong Zhang  */
124b561aa9dSHong Zhang #undef __FUNCT__
12525023028SHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable"
126bf9555e6SHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
127e9e4536cSHong Zhang {
128e9e4536cSHong Zhang   PetscErrorCode     ierr;
129bf9555e6SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
130bf9555e6SHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data;
131bf9555e6SHong Zhang   const PetscInt     *Ai=a->i,*Aj=a->j,*Bi=b->i,*Bj=b->j,am=A->rmap->N,bm=B->rmap->N,bn=B->cmap->N;
132bf9555e6SHong Zhang   const PetscInt     *aj=Aj,*bi=Bi,*bj;
133971236abSHong Zhang   PetscInt           *ci,*cj;
134bf9555e6SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,crmax,ndouble=0;
135e9e4536cSHong Zhang   PetscBT            bt;
136bf9555e6SHong Zhang   PetscInt           nlnk_max,lnk_max=bn,*lnk,*nlnk;
137e9e4536cSHong Zhang 
138e9e4536cSHong Zhang   PetscFunctionBegin;
139bf9555e6SHong Zhang   /* Allocate ci array, arrays for fill computation and */
140bf9555e6SHong Zhang   /* free space for accumulating nonzero column info */
141e9e4536cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
142e9e4536cSHong Zhang   ci[0] = 0;
143e9e4536cSHong Zhang 
144bf9555e6SHong Zhang   /* create and initialize a condensed linked list */
145bf9555e6SHong Zhang   crmax = a->rmax*b->rmax;
146bf9555e6SHong Zhang   nlnk_max = PetscMin(crmax,bn);
147164eb823SHong Zhang #if defined(DEBUG_MATMATMULT)
148bf9555e6SHong Zhang   printf("LLCondensedCreate nlnk_max=%d, bn %d, crmax %d\n",nlnk_max,bn,crmax);
149164eb823SHong Zhang #endif
150bf9555e6SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,lnk_max,lnk,nlnk,bt);CHKERRQ(ierr);
151bf9555e6SHong Zhang 
152bf9555e6SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
153bf9555e6SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(Ai[am]+Bi[bm])),&free_space);CHKERRQ(ierr);
154bf9555e6SHong Zhang   current_space = free_space;
155bf9555e6SHong Zhang 
156bf9555e6SHong Zhang   /* Determine ci and cj */
157e9e4536cSHong Zhang   for (i=0; i<am; i++) {
15831e66dd0SHong Zhang     anzi = Ai[i+1] - Ai[i];
159e9e4536cSHong Zhang     cnzi = 0;
16031e66dd0SHong Zhang     aj   = Aj + Ai[i];
161e9e4536cSHong Zhang     for (j=0; j<anzi; j++){
162e9e4536cSHong Zhang       brow = aj[j];
163e9e4536cSHong Zhang       bnzj = bi[brow+1] - bi[brow];
16431e66dd0SHong Zhang       bj   = Bj + bi[brow];
165bf9555e6SHong Zhang       /* add non-zero cols of B into the condensed sorted linked list lnk */
166bf9555e6SHong Zhang       ierr = PetscLLCondensedAddSorted(nlnk_max,lnk_max,bnzj,bj,nlnk,lnk,bt);CHKERRQ(ierr);
167e9e4536cSHong Zhang     }
1688d298b4fSHong Zhang     cnzi    = *nlnk;
169e9e4536cSHong Zhang     ci[i+1] = ci[i] + cnzi;
170bf9555e6SHong Zhang 
171bf9555e6SHong Zhang     /* If free space is not available, make more free space */
172bf9555e6SHong Zhang     /* Double the amount of total space in the list */
173bf9555e6SHong Zhang     if (current_space->local_remaining<cnzi) {
174bf9555e6SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
175bf9555e6SHong Zhang       ndouble++;
176e9e4536cSHong Zhang     }
177bf9555e6SHong Zhang 
178bf9555e6SHong Zhang     /* Copy data into free space, then initialize lnk */
179bf9555e6SHong Zhang     ierr = PetscLLCondensedClean(nlnk_max,lnk_max,cnzi,current_space->array,nlnk,lnk,bt);CHKERRQ(ierr);
180bf9555e6SHong Zhang     current_space->array           += cnzi;
181bf9555e6SHong Zhang     current_space->local_used      += cnzi;
182bf9555e6SHong Zhang     current_space->local_remaining -= cnzi;
183bf9555e6SHong Zhang     ci[i+1] = ci[i] + cnzi;
184bf9555e6SHong Zhang   }
185bf9555e6SHong Zhang 
186bf9555e6SHong Zhang   /* Column indices are in the list of free space */
1878d298b4fSHong Zhang   /* Allocate and copy column indices to cj, then destroy list of free space, lnk and bt */
188e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
189bf9555e6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
190bf9555e6SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,bt);CHKERRQ(ierr);
191e9e4536cSHong Zhang 
192e9e4536cSHong Zhang   *Ci           = ci;
193e9e4536cSHong Zhang   *Cj           = cj;
194bf9555e6SHong Zhang   *nspacedouble = ndouble;
195e9e4536cSHong Zhang   PetscFunctionReturn(0);
196e9e4536cSHong Zhang }
197e9e4536cSHong Zhang 
198e9e4536cSHong Zhang #undef __FUNCT__
199b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
200b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
201b561aa9dSHong Zhang {
202b561aa9dSHong Zhang   PetscErrorCode ierr;
203b561aa9dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
204971236abSHong Zhang   PetscInt       *ai=a->i,*bi=b->i,*ci,*cj;
205b561aa9dSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
206b561aa9dSHong Zhang   MatScalar      *ca;
207b561aa9dSHong Zhang   PetscReal      afill;
208b561aa9dSHong Zhang 
209b561aa9dSHong Zhang   PetscFunctionBegin;
210e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
211971236abSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultSymbolic_SeqAIJ_SeqAIJ...\n");CHKERRQ(ierr);
212e9e4536cSHong Zhang #endif
213971236abSHong Zhang   /* Get ci and cj */
214971236abSHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(A,B,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
215b561aa9dSHong Zhang 
21658c24d83SHong Zhang   /* Allocate space for ca */
21758c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
21858c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
21958c24d83SHong Zhang 
22026be0446SHong Zhang   /* put together the new symbolic matrix */
2217adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
22258c24d83SHong Zhang 
22358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
22458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
22558c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
226e6b907acSBarry Smith   c->free_a   = PETSC_TRUE;
227e6b907acSBarry Smith   c->free_ij  = PETSC_TRUE;
22858c24d83SHong Zhang   c->nonew    = 0;
2298cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ;
23058c24d83SHong Zhang 
23136ec6d2dSHong Zhang   ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr);
23225023028SHong Zhang   ierr = PetscMemzero(c->matmult_abdense,bn*sizeof(PetscScalar));CHKERRQ(ierr);
233c58ca2e3SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */
2340b7e3e3dSHong Zhang 
2357212da7cSHong Zhang   /* set MatInfo */
2367212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
237f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
2387212da7cSHong Zhang   c->maxnz                     = ci[am];
2397212da7cSHong Zhang   c->nz                        = ci[am];
2407212da7cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
2417212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
2427212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
2437212da7cSHong Zhang 
2447212da7cSHong Zhang #if defined(PETSC_USE_INFO)
2457212da7cSHong Zhang   if (ci[am]) {
246f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
247f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
248f2b054eeSHong Zhang   } else {
249f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
250be0fcf8dSHong Zhang   }
251f2b054eeSHong Zhang #endif
25258c24d83SHong Zhang   PetscFunctionReturn(0);
25358c24d83SHong Zhang }
254d50806bdSBarry Smith 
25526be0446SHong Zhang #undef __FUNCT__
25626be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
257dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
258d50806bdSBarry Smith {
259dfbe8321SBarry Smith   PetscErrorCode ierr;
260d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
261d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
262d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
263d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
26438baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
265c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
266519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
26736ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
2680b7e3e3dSHong Zhang   PetscScalar    *ab_dense=c->matmult_abdense;
269d50806bdSBarry Smith 
270d50806bdSBarry Smith   PetscFunctionBegin;
271c124e916SHong Zhang   /* clean old values in C */
272c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
273d50806bdSBarry Smith   /* Traverse A row-wise. */
274d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
275d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
276d50806bdSBarry Smith   for (i=0; i<am; i++) {
277d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
278d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
279519eb980SHong Zhang       brow = aj[j];
280d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
281d50806bdSBarry Smith       bjj  = bj + bi[brow];
282d50806bdSBarry Smith       baj  = ba + bi[brow];
283519eb980SHong Zhang       /* perform dense axpy */
28436ec6d2dSHong Zhang       valtmp = aa[j];
285519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
28636ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
287519eb980SHong Zhang       }
288519eb980SHong Zhang       flops += 2*bnzi;
289519eb980SHong Zhang     }
290c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
291c58ca2e3SHong Zhang 
292c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
293519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
294519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
295519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
296519eb980SHong Zhang     }
297519eb980SHong Zhang     flops += cnzi;
298519eb980SHong Zhang     cj += cnzi; ca += cnzi;
299519eb980SHong Zhang   }
300c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
301c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
302c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
303c58ca2e3SHong Zhang   PetscFunctionReturn(0);
304c58ca2e3SHong Zhang }
305c58ca2e3SHong Zhang 
306c58ca2e3SHong Zhang #undef __FUNCT__
30725023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
30825023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
309c58ca2e3SHong Zhang {
310c58ca2e3SHong Zhang   PetscErrorCode ierr;
311c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
312c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
313c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
314c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
315c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
316c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
317c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
31836ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
319c58ca2e3SHong Zhang   PetscInt       nextb;
320c58ca2e3SHong Zhang 
321c58ca2e3SHong Zhang   PetscFunctionBegin;
322e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
323971236abSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr);
324e9e4536cSHong Zhang #endif
325c58ca2e3SHong Zhang   /* clean old values in C */
326c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
327c58ca2e3SHong Zhang   /* Traverse A row-wise. */
328c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
329c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
330519eb980SHong Zhang   for (i=0;i<am;i++) {
331519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
332519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
333519eb980SHong Zhang     for (j=0;j<anzi;j++) {
334519eb980SHong Zhang       brow = aj[j];
335519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
336519eb980SHong Zhang       bjj  = bj + bi[brow];
337519eb980SHong Zhang       baj  = ba + bi[brow];
338519eb980SHong Zhang       /* perform sparse axpy */
33936ec6d2dSHong Zhang       valtmp = aa[j];
340c124e916SHong Zhang       nextb  = 0;
341c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
342c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
34336ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
344c124e916SHong Zhang         }
345d50806bdSBarry Smith       }
346d50806bdSBarry Smith       flops += 2*bnzi;
347d50806bdSBarry Smith     }
348519eb980SHong Zhang     aj += anzi; aa += anzi;
349519eb980SHong Zhang     cj += cnzi; ca += cnzi;
350d50806bdSBarry Smith   }
351c58ca2e3SHong Zhang 
352716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
353716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
354d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
355d50806bdSBarry Smith   PetscFunctionReturn(0);
356d50806bdSBarry Smith }
357bc011b1eSHong Zhang 
358e9e4536cSHong Zhang #undef __FUNCT__
35925023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
36025023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
361e9e4536cSHong Zhang {
362e9e4536cSHong Zhang   PetscErrorCode ierr;
363e9e4536cSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
364bf9555e6SHong Zhang   PetscInt       *ai=a->i,*bi=b->i,*ci,*cj;
365e9e4536cSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
366e9e4536cSHong Zhang   MatScalar      *ca;
367e9e4536cSHong Zhang   PetscReal      afill;
368e9e4536cSHong Zhang 
369e9e4536cSHong Zhang   PetscFunctionBegin;
370e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
371164eb823SHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable Armax %d, Brmax %d\n",a->rmax,b->rmax);CHKERRQ(ierr);
372e9e4536cSHong Zhang #endif
373e9e4536cSHong Zhang   /* Get ci and cj */
374bf9555e6SHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable(A,B,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
375e9e4536cSHong Zhang 
376e9e4536cSHong Zhang   /* Allocate space for ca */
377e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
378e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
379e9e4536cSHong Zhang 
380e9e4536cSHong Zhang   /* put together the new symbolic matrix */
381e9e4536cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
382e9e4536cSHong Zhang 
383e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
384e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
385e9e4536cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
386e9e4536cSHong Zhang   c->free_a   = PETSC_TRUE;
387e9e4536cSHong Zhang   c->free_ij  = PETSC_TRUE;
388e9e4536cSHong Zhang   c->nonew    = 0;
38925023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
390e9e4536cSHong Zhang 
391e9e4536cSHong Zhang   /* set MatInfo */
392e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
393e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
394e9e4536cSHong Zhang   c->maxnz                     = ci[am];
395e9e4536cSHong Zhang   c->nz                        = ci[am];
396e9e4536cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
397e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
398e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
399e9e4536cSHong Zhang 
400e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
401e9e4536cSHong Zhang   if (ci[am]) {
402e9e4536cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
403e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
404e9e4536cSHong Zhang   } else {
405e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
406e9e4536cSHong Zhang   }
407e9e4536cSHong Zhang #endif
408e9e4536cSHong Zhang   PetscFunctionReturn(0);
409e9e4536cSHong Zhang }
410e9e4536cSHong Zhang 
411e9e4536cSHong Zhang 
412d2853540SHong Zhang /* This routine is not used. Should be removed! */
413bc011b1eSHong Zhang #undef __FUNCT__
4146fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
4156fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4165df89d91SHong Zhang {
417bc011b1eSHong Zhang   PetscErrorCode ierr;
418bc011b1eSHong Zhang 
419bc011b1eSHong Zhang   PetscFunctionBegin;
420bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
4216fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
422bc011b1eSHong Zhang   }
4236fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
424bc011b1eSHong Zhang   PetscFunctionReturn(0);
425bc011b1eSHong Zhang }
426bc011b1eSHong Zhang 
427bc011b1eSHong Zhang #undef __FUNCT__
428e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
429e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
4302128a86cSHong Zhang {
4312128a86cSHong Zhang   PetscErrorCode      ierr;
432e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
4332128a86cSHong Zhang 
4342128a86cSHong Zhang   PetscFunctionBegin;
435b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
4362128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
4372128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
4382128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
4392128a86cSHong Zhang   PetscFunctionReturn(0);
4402128a86cSHong Zhang }
4412128a86cSHong Zhang 
4422128a86cSHong Zhang #undef __FUNCT__
4432128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
4442128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
4452128a86cSHong Zhang {
4462128a86cSHong Zhang   PetscErrorCode      ierr;
4472128a86cSHong Zhang   PetscContainer      container;
448e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
4492128a86cSHong Zhang 
4502128a86cSHong Zhang   PetscFunctionBegin;
451e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
4522128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
4532128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
4542128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
4552128a86cSHong Zhang   if (A->ops->destroy) {
4562128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
4572128a86cSHong Zhang   }
458e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
4592128a86cSHong Zhang   PetscFunctionReturn(0);
4602128a86cSHong Zhang }
4612128a86cSHong Zhang 
4622128a86cSHong Zhang #undef __FUNCT__
4636fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
4646fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
465bc011b1eSHong Zhang {
4665df89d91SHong Zhang   PetscErrorCode      ierr;
46781d82fe4SHong Zhang   Mat                 Bt;
46881d82fe4SHong Zhang   PetscInt            *bti,*btj;
469e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
4702128a86cSHong Zhang   PetscContainer      container;
471d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
472d2853540SHong Zhang 
47381d82fe4SHong Zhang   PetscFunctionBegin;
474d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
47581d82fe4SHong Zhang    /* create symbolic Bt */
47681d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
47781d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
47881d82fe4SHong Zhang 
47981d82fe4SHong Zhang   /* get symbolic C=A*Bt */
48081d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
48181d82fe4SHong Zhang 
4822128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
483e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
4842128a86cSHong Zhang 
4852128a86cSHong Zhang   /* attach the supporting struct to C */
4862128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
4872128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
488e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
489e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
4902128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4912128a86cSHong Zhang 
4922128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
4932128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
4942128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
4952128a86cSHong Zhang 
496d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
497d2853540SHong Zhang   etime2 += tf - t0;
498d2853540SHong Zhang 
499f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
5002128a86cSHong Zhang   if (multtrans->usecoloring){
501b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
502b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
5032128a86cSHong Zhang     ISColoring           iscoloring;
5042128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
505d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
506e8354b3eSHong Zhang 
507e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
508502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
509d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
510d2853540SHong Zhang     etime0 += tf - t0;
511d2853540SHong Zhang 
512d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
513b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
5142128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
5152128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
516e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
517d2853540SHong Zhang     etime01 += tf - t0;
5182128a86cSHong Zhang 
519e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
5202128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
5212128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
5222128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
5232128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
5242128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
5252128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
5262128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
5272128a86cSHong Zhang 
5282128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
5292128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
5302128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
5312128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
5322128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
5332128a86cSHong Zhang     multtrans->ABt_den = C_dense;
534e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
535e8354b3eSHong Zhang     etime1 += tf - t0;
536f94ccd6cSHong Zhang 
537f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
5381ce71dffSSatish Balay     {
539f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
540f94ccd6cSHong 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));
541f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
5421ce71dffSSatish Balay     }
543f94ccd6cSHong Zhang #endif
5442128a86cSHong Zhang   }
545*e99be685SHong Zhang   /* clean up */
546*e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
547*e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
5482128a86cSHong Zhang 
549f94ccd6cSHong Zhang 
550f94ccd6cSHong Zhang 
55181d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
55281d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
5531fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
5541fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5551fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
5561fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
5571fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
5581fa1209cSHong Zhang   MatScalar          *ca;
5591fa1209cSHong Zhang   PetscBT            lnkbt;
5601fa1209cSHong Zhang   PetscReal          afill;
5615df89d91SHong Zhang 
5621fa1209cSHong Zhang   /* Allocate row pointer array ci  */
5631fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5641fa1209cSHong Zhang   ci[0] = 0;
5651fa1209cSHong Zhang 
5661fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
5671fa1209cSHong Zhang   nlnk = bm+1;
5681fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5691fa1209cSHong Zhang 
5701fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
5711fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5721fa1209cSHong Zhang   current_space = free_space;
5731fa1209cSHong Zhang 
5741fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
5751fa1209cSHong Zhang   for (i=0; i<am; i++) {
5761fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
5771fa1209cSHong Zhang     cnzi = 0;
5781fa1209cSHong Zhang     acol = aj + ai[i];
5791fa1209cSHong Zhang     for (j=0; j<bm; j++){
5801fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
5811fa1209cSHong Zhang       bcol= bj + bi[j];
58281d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
5831fa1209cSHong Zhang       ka = 0; kb = 0;
5841fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
58581d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
58681d82fe4SHong Zhang         if (ka == anzi) break;
58781d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
58881d82fe4SHong Zhang         if (kb == bnzj) break;
5891fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
5901fa1209cSHong Zhang           index[0] = j;
5911fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5921fa1209cSHong Zhang           cnzi++;
5931fa1209cSHong Zhang           break;
5941fa1209cSHong Zhang         }
5951fa1209cSHong Zhang       }
5961fa1209cSHong Zhang     }
5971fa1209cSHong Zhang 
5981fa1209cSHong Zhang     /* If free space is not available, make more free space */
5991fa1209cSHong Zhang     /* Double the amount of total space in the list */
6001fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
6011fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
6021fa1209cSHong Zhang       nspacedouble++;
6031fa1209cSHong Zhang     }
6041fa1209cSHong Zhang 
6051fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
6061fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
6071fa1209cSHong Zhang     current_space->array           += cnzi;
6081fa1209cSHong Zhang     current_space->local_used      += cnzi;
6091fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
6101fa1209cSHong Zhang 
6111fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
6121fa1209cSHong Zhang   }
6131fa1209cSHong Zhang 
6141fa1209cSHong Zhang 
6151fa1209cSHong Zhang   /* Column indices are in the list of free space.
6161fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
6171fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6181fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6191fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
6201fa1209cSHong Zhang 
6211fa1209cSHong Zhang   /* Allocate space for ca */
6221fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
6231fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
6241fa1209cSHong Zhang 
6251fa1209cSHong Zhang   /* put together the new symbolic matrix */
6261fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
6271fa1209cSHong Zhang 
6281fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6291fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6301fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
6311fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
6321fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
6331fa1209cSHong Zhang   c->nonew    = 0;
6341fa1209cSHong Zhang 
6351fa1209cSHong Zhang   /* set MatInfo */
6361fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6371fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
6381fa1209cSHong Zhang   c->maxnz                     = ci[am];
6391fa1209cSHong Zhang   c->nz                        = ci[am];
6401fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
6411fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
6421fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
6431fa1209cSHong Zhang 
6441fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
6451fa1209cSHong Zhang   if (ci[am]) {
6461fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
6476fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
6481fa1209cSHong Zhang   } else {
6491fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6501fa1209cSHong Zhang   }
6511fa1209cSHong Zhang #endif
65281d82fe4SHong Zhang #endif
6535df89d91SHong Zhang   PetscFunctionReturn(0);
6545df89d91SHong Zhang }
6555df89d91SHong Zhang 
6566973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
6575df89d91SHong Zhang #undef __FUNCT__
6586fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
6596fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
6605df89d91SHong Zhang {
6615df89d91SHong Zhang   PetscErrorCode ierr;
6625df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
663e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
66481d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
6655df89d91SHong Zhang   PetscLogDouble flops=0.0;
6666973769bSHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
667e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
6682128a86cSHong Zhang   PetscContainer      container;
6696973769bSHong Zhang #if defined(USE_ARRAY)
6706973769bSHong Zhang   MatScalar      *spdot;
6716973769bSHong Zhang #endif
6725df89d91SHong Zhang 
6735df89d91SHong Zhang   PetscFunctionBegin;
674e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
6752128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
6762128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
6772128a86cSHong Zhang   if (multtrans->usecoloring){
678b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
679c8db22f6SHong Zhang     Mat                   Bt_dense;
680c8db22f6SHong Zhang     PetscInt              m,n;
681b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
682a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
683a0b95ffbSSatish Balay 
6842128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
685c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
686c8db22f6SHong Zhang 
687b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
6882128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
689b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
6902128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
691b2d2b04fSHong Zhang     etime0 += tf - t0;
692c8db22f6SHong Zhang 
693c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
6942128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
6952128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
6962128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
6972128a86cSHong Zhang     etime2 += tf - t0;
698c8db22f6SHong Zhang 
6992128a86cSHong Zhang     /* Recover C from C_dense */
7002128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
701b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
7022128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7032128a86cSHong Zhang     etime1 += tf - t0;
704f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
705f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
706f94ccd6cSHong Zhang #endif
707c8db22f6SHong Zhang     PetscFunctionReturn(0);
708c8db22f6SHong Zhang   }
709c8db22f6SHong Zhang 
7106973769bSHong Zhang #if defined(USE_ARRAY)
7116973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
712e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
713e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
7146973769bSHong Zhang #endif
7156973769bSHong Zhang 
71681d82fe4SHong Zhang   /* clear old values in C */
71781d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
7185df89d91SHong Zhang 
7191fa1209cSHong Zhang   for (i=0; i<cm; i++) {
72081d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
72181d82fe4SHong Zhang     acol = aj + ai[i];
7226973769bSHong Zhang     aval = aa + ai[i];
7231fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
7241fa1209cSHong Zhang     ccol = cj + ci[i];
7256973769bSHong Zhang     cval = ca + ci[i];
7261fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
72781d82fe4SHong Zhang       brow = ccol[j];
72881d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
72981d82fe4SHong Zhang       bcol = bj + bi[brow];
7306973769bSHong Zhang       bval = ba + bi[brow];
7316973769bSHong Zhang 
73281d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
7336973769bSHong Zhang #if defined(USE_ARRAY)
7346973769bSHong Zhang       /* put ba to spdot array */
7356973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
7366973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
7376973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
7386973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
7396973769bSHong Zhang       }
7406973769bSHong Zhang       /* zero spdot array */
7416973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
7426973769bSHong Zhang #else
74381d82fe4SHong Zhang       nexta = 0; nextb = 0;
74481d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
74581d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
74681d82fe4SHong Zhang         if (nexta == anzi) break;
74781d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
74881d82fe4SHong Zhang         if (nextb == bnzj) break;
74981d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
7506973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
75181d82fe4SHong Zhang           nexta++; nextb++;
75281d82fe4SHong Zhang           flops += 2;
7531fa1209cSHong Zhang         }
7541fa1209cSHong Zhang       }
7556973769bSHong Zhang #endif
75681d82fe4SHong Zhang     }
75781d82fe4SHong Zhang   }
75881d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
75981d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76081d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
7616973769bSHong Zhang #if defined(USE_ARRAY)
7626973769bSHong Zhang   ierr = PetscFree(spdot);
7636973769bSHong Zhang #endif
7645df89d91SHong Zhang   PetscFunctionReturn(0);
7655df89d91SHong Zhang }
7665df89d91SHong Zhang 
7675df89d91SHong Zhang #undef __FUNCT__
76875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
76975648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
7705df89d91SHong Zhang   PetscErrorCode ierr;
7715df89d91SHong Zhang 
7725df89d91SHong Zhang   PetscFunctionBegin;
7735df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
77475648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
7755df89d91SHong Zhang   }
77675648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
7775df89d91SHong Zhang   PetscFunctionReturn(0);
7785df89d91SHong Zhang }
7795df89d91SHong Zhang 
7805df89d91SHong Zhang #undef __FUNCT__
78175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
78275648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
7835df89d91SHong Zhang {
784bc011b1eSHong Zhang   PetscErrorCode ierr;
785bc011b1eSHong Zhang   Mat            At;
78638baddfdSBarry Smith   PetscInt       *ati,*atj;
787bc011b1eSHong Zhang 
788bc011b1eSHong Zhang   PetscFunctionBegin;
789bc011b1eSHong Zhang   /* create symbolic At */
790bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
791d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
792bc011b1eSHong Zhang 
793bc011b1eSHong Zhang   /* get symbolic C=At*B */
794bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
795bc011b1eSHong Zhang 
796bc011b1eSHong Zhang   /* clean up */
7976bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
798bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
799bc011b1eSHong Zhang   PetscFunctionReturn(0);
800bc011b1eSHong Zhang }
801bc011b1eSHong Zhang 
802bc011b1eSHong Zhang #undef __FUNCT__
80375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
80475648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
805bc011b1eSHong Zhang {
806bc011b1eSHong Zhang   PetscErrorCode ierr;
8070fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
808d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
809d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
810d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
8110fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
812bc011b1eSHong Zhang 
813bc011b1eSHong Zhang   PetscFunctionBegin;
814bc011b1eSHong Zhang   /* clear old values in C */
815bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
816bc011b1eSHong Zhang 
817bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
818bc011b1eSHong Zhang   for (i=0;i<am;i++) {
819bc011b1eSHong Zhang     bj   = b->j + bi[i];
820bc011b1eSHong Zhang     ba   = b->a + bi[i];
821bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
822bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
823bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
824bc011b1eSHong Zhang       nextb = 0;
8250fbc74f4SHong Zhang       crow  = *aj++;
826bc011b1eSHong Zhang       cjj   = cj + ci[crow];
827bc011b1eSHong Zhang       caj   = ca + ci[crow];
828bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
829bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
8300fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
8310fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
832bc011b1eSHong Zhang           nextb++;
833bc011b1eSHong Zhang         }
834bc011b1eSHong Zhang       }
835bc011b1eSHong Zhang       flops += 2*bnzi;
8360fbc74f4SHong Zhang       aa++;
837bc011b1eSHong Zhang     }
838bc011b1eSHong Zhang   }
839bc011b1eSHong Zhang 
840bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
841bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
842bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
843bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
844bc011b1eSHong Zhang   PetscFunctionReturn(0);
845bc011b1eSHong Zhang }
8469123193aSHong Zhang 
847ec01deb9SMatthew Knepley EXTERN_C_BEGIN
8489123193aSHong Zhang #undef __FUNCT__
8499123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
8509123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8519123193aSHong Zhang {
8529123193aSHong Zhang   PetscErrorCode ierr;
8539123193aSHong Zhang 
8549123193aSHong Zhang   PetscFunctionBegin;
8559123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
8569123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
8579123193aSHong Zhang   }
8589123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
8599123193aSHong Zhang   PetscFunctionReturn(0);
8609123193aSHong Zhang }
861ec01deb9SMatthew Knepley EXTERN_C_END
862edd81eecSMatthew Knepley 
8639123193aSHong Zhang #undef __FUNCT__
8649123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
8659123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
8669123193aSHong Zhang {
8679123193aSHong Zhang   PetscErrorCode ierr;
8689123193aSHong Zhang 
8699123193aSHong Zhang   PetscFunctionBegin;
8705a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
8718cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
8729123193aSHong Zhang   PetscFunctionReturn(0);
8739123193aSHong Zhang }
8749123193aSHong Zhang 
8759123193aSHong Zhang #undef __FUNCT__
8769123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
8779123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
8789123193aSHong Zhang {
879f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
8809123193aSHong Zhang   PetscErrorCode ierr;
881dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
882dd6ea824SBarry Smith   MatScalar      *aa;
883d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
884f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
8859123193aSHong Zhang 
8869123193aSHong Zhang   PetscFunctionBegin;
887f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
888e32f2f54SBarry 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);
889e32f2f54SBarry 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);
890e32f2f54SBarry 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);
891f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
892f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
893f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
894f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
895f32d5d43SBarry Smith     colam = col*am;
896f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
897f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
898f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
899f32d5d43SBarry Smith       aj  = a->j + a->i[i];
900f32d5d43SBarry Smith       aa  = a->a + a->i[i];
901f32d5d43SBarry Smith       for (j=0; j<n; j++) {
902f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
903f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
904f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
905f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
9069123193aSHong Zhang       }
907f32d5d43SBarry Smith       c[colam + i]       = r1;
908f32d5d43SBarry Smith       c[colam + am + i]  = r2;
909f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
910f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
911f32d5d43SBarry Smith     }
912f32d5d43SBarry Smith     b1 += bm4;
913f32d5d43SBarry Smith     b2 += bm4;
914f32d5d43SBarry Smith     b3 += bm4;
915f32d5d43SBarry Smith     b4 += bm4;
916f32d5d43SBarry Smith   }
917f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
918f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
919f32d5d43SBarry Smith       r1 = 0.0;
920f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
921f32d5d43SBarry Smith       aj  = a->j + a->i[i];
922f32d5d43SBarry Smith       aa  = a->a + a->i[i];
923f32d5d43SBarry Smith 
924f32d5d43SBarry Smith       for (j=0; j<n; j++) {
925f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
926f32d5d43SBarry Smith       }
927f32d5d43SBarry Smith       c[col*am + i]     = r1;
928f32d5d43SBarry Smith     }
929f32d5d43SBarry Smith     b1 += bm;
930f32d5d43SBarry Smith   }
931dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
932f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
933f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
934f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
935f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
936f32d5d43SBarry Smith   PetscFunctionReturn(0);
937f32d5d43SBarry Smith }
938f32d5d43SBarry Smith 
939f32d5d43SBarry Smith /*
9404324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
941f32d5d43SBarry Smith */
942f32d5d43SBarry Smith #undef __FUNCT__
943f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
944f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
945f32d5d43SBarry Smith {
946f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
947f32d5d43SBarry Smith   PetscErrorCode ierr;
948dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
949dd6ea824SBarry Smith   MatScalar      *aa;
950d0f46423SBarry 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;
9514324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
952f32d5d43SBarry Smith 
953f32d5d43SBarry Smith   PetscFunctionBegin;
954f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
955f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
956f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
957f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
9584324174eSBarry Smith 
9594324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
9604324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
9614324174eSBarry Smith       colam = col*am;
9624324174eSBarry Smith       arm   = a->compressedrow.nrows;
9634324174eSBarry Smith       ii    = a->compressedrow.i;
9644324174eSBarry Smith       ridx  = a->compressedrow.rindex;
9654324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
9664324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
9674324174eSBarry Smith 	n   = ii[i+1] - ii[i];
9684324174eSBarry Smith 	aj  = a->j + ii[i];
9694324174eSBarry Smith 	aa  = a->a + ii[i];
9704324174eSBarry Smith 	for (j=0; j<n; j++) {
9714324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
9724324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
9734324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
9744324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
9754324174eSBarry Smith 	}
9764324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
9774324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
9784324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
9794324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
9804324174eSBarry Smith       }
9814324174eSBarry Smith       b1 += bm4;
9824324174eSBarry Smith       b2 += bm4;
9834324174eSBarry Smith       b3 += bm4;
9844324174eSBarry Smith       b4 += bm4;
9854324174eSBarry Smith     }
9864324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
9874324174eSBarry Smith       colam = col*am;
9884324174eSBarry Smith       arm   = a->compressedrow.nrows;
9894324174eSBarry Smith       ii    = a->compressedrow.i;
9904324174eSBarry Smith       ridx  = a->compressedrow.rindex;
9914324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
9924324174eSBarry Smith 	r1 = 0.0;
9934324174eSBarry Smith 	n   = ii[i+1] - ii[i];
9944324174eSBarry Smith 	aj  = a->j + ii[i];
9954324174eSBarry Smith 	aa  = a->a + ii[i];
9964324174eSBarry Smith 
9974324174eSBarry Smith 	for (j=0; j<n; j++) {
9984324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
9994324174eSBarry Smith 	}
10004324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
10014324174eSBarry Smith       }
10024324174eSBarry Smith       b1 += bm;
10034324174eSBarry Smith     }
10044324174eSBarry Smith   } else {
1005f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
1006f32d5d43SBarry Smith       colam = col*am;
1007f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1008f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
1009f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1010f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1011f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1012f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1013f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
1014f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
1015f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
1016f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
1017f32d5d43SBarry Smith 	}
1018f32d5d43SBarry Smith 	c[colam + i]       += r1;
1019f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
1020f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
1021f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
1022f32d5d43SBarry Smith       }
1023f32d5d43SBarry Smith       b1 += bm4;
1024f32d5d43SBarry Smith       b2 += bm4;
1025f32d5d43SBarry Smith       b3 += bm4;
1026f32d5d43SBarry Smith       b4 += bm4;
1027f32d5d43SBarry Smith     }
1028f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
1029f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1030f32d5d43SBarry Smith 	r1 = 0.0;
1031f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1032f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1033f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1034f32d5d43SBarry Smith 
1035f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1036f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
1037f32d5d43SBarry Smith 	}
1038f32d5d43SBarry Smith 	c[col*am + i]     += r1;
1039f32d5d43SBarry Smith       }
1040f32d5d43SBarry Smith       b1 += bm;
1041f32d5d43SBarry Smith     }
10424324174eSBarry Smith   }
1043dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
1044f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1045f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
10469123193aSHong Zhang   PetscFunctionReturn(0);
10479123193aSHong Zhang }
1048b1683b59SHong Zhang 
1049b1683b59SHong Zhang #undef __FUNCT__
1050b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1051b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1052c8db22f6SHong Zhang {
1053c8db22f6SHong Zhang   PetscErrorCode ierr;
10542f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
10552f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
10562f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
10572f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
10582f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
10592f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1060c8db22f6SHong Zhang 
1061c8db22f6SHong Zhang   PetscFunctionBegin;
10622f5f1f90SHong Zhang   btval_den=btdense->v;
10632f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
10642f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
10652f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
10662f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1067d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
10682f5f1f90SHong Zhang       btcol = bj + bi[col];
10692f5f1f90SHong Zhang       btval = ba + bi[col];
10702f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1071d2853540SHong Zhang       for (j=0; j<anz; j++){
10722f5f1f90SHong Zhang         brow            = btcol[j];
10732f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1074c8db22f6SHong Zhang       }
1075c8db22f6SHong Zhang     }
10762f5f1f90SHong Zhang     btval_den += m;
1077c8db22f6SHong Zhang   }
1078c8db22f6SHong Zhang   PetscFunctionReturn(0);
1079c8db22f6SHong Zhang }
1080c8db22f6SHong Zhang 
1081c8db22f6SHong Zhang #undef __FUNCT__
1082b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1083b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
10848972f759SHong Zhang {
10858972f759SHong Zhang   PetscErrorCode ierr;
1086b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
10872f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1088b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
10892f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
10908972f759SHong Zhang 
10918972f759SHong Zhang   PetscFunctionBegin;
10928972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
1093b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
1094b2d2b04fSHong Zhang   cp_den = ca_den;
10952f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
10962f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
10972f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
10982f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
10992f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
11002f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1101b2d2b04fSHong Zhang     }
1102b2d2b04fSHong Zhang     cp_den += m;
1103b2d2b04fSHong Zhang   }
1104b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
11058972f759SHong Zhang   PetscFunctionReturn(0);
11068972f759SHong Zhang }
11078972f759SHong Zhang 
1108*e99be685SHong Zhang /*
1109*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1110*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1111*e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1112*e99be685SHong Zhang  */
1113*e99be685SHong Zhang #undef __FUNCT__
1114*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
1115*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1116*e99be685SHong Zhang {
1117*e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1118*e99be685SHong Zhang   PetscErrorCode ierr;
1119*e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1120*e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1121*e99be685SHong Zhang   PetscInt       *cspidx;
1122*e99be685SHong Zhang 
1123*e99be685SHong Zhang   PetscFunctionBegin;
1124*e99be685SHong Zhang   *nn = n;
1125*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1126*e99be685SHong Zhang   if (symmetric) {
1127*e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
1128*e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
1129*e99be685SHong Zhang   } else {
1130*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1131*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1132*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1133*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1134*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1135*e99be685SHong Zhang     jj = a->j;
1136*e99be685SHong Zhang     for (i=0; i<nz; i++) {
1137*e99be685SHong Zhang       collengths[jj[i]]++;
1138*e99be685SHong Zhang     }
1139*e99be685SHong Zhang     cia[0] = oshift;
1140*e99be685SHong Zhang     for (i=0; i<n; i++) {
1141*e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1142*e99be685SHong Zhang     }
1143*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1144*e99be685SHong Zhang     jj   = a->j;
1145*e99be685SHong Zhang     for (row=0; row<m; row++) {
1146*e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1147*e99be685SHong Zhang       for (i=0; i<mr; i++) {
1148*e99be685SHong Zhang         col = *jj++;
1149*e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1150*e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1151*e99be685SHong Zhang       }
1152*e99be685SHong Zhang     }
1153*e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1154*e99be685SHong Zhang     *ia = cia; *ja = cja;
1155*e99be685SHong Zhang     *spidx = cspidx;
1156*e99be685SHong Zhang   }
1157*e99be685SHong Zhang   PetscFunctionReturn(0);
1158*e99be685SHong Zhang }
1159*e99be685SHong Zhang 
1160*e99be685SHong Zhang #undef __FUNCT__
1161*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1162*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1163*e99be685SHong Zhang {
1164*e99be685SHong Zhang   PetscErrorCode ierr;
1165*e99be685SHong Zhang 
1166*e99be685SHong Zhang   PetscFunctionBegin;
1167*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1168*e99be685SHong Zhang 
1169*e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1170*e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1171*e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1172*e99be685SHong Zhang   PetscFunctionReturn(0);
1173*e99be685SHong Zhang }
1174*e99be685SHong Zhang 
11758972f759SHong Zhang #undef __FUNCT__
1176b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1177b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1178b1683b59SHong Zhang {
1179b1683b59SHong Zhang   PetscErrorCode ierr;
1180d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1181b1683b59SHong Zhang   const PetscInt *is;
1182b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1183b1683b59SHong Zhang   IS             *isa;
1184b1683b59SHong Zhang   PetscBool      done;
1185b1683b59SHong Zhang   PetscBool      flg1,flg2;
1186b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1187*e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1188d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1189b1683b59SHong Zhang 
1190b1683b59SHong Zhang   PetscFunctionBegin;
1191b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1192*e99be685SHong Zhang 
1193b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1194b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1195b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1196b1683b59SHong Zhang   if (flg1 || flg2) {
1197b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1198b1683b59SHong Zhang   }
1199b1683b59SHong Zhang 
1200b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1201b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1202b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1203b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1204b1683b59SHong Zhang   c->rstart = 0;
1205b1683b59SHong Zhang 
1206b1683b59SHong Zhang   c->ncolors = nis;
1207b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1208b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1209d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1210d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1211d2853540SHong Zhang   colorforrow[0]    = 0;
1212d2853540SHong Zhang   rows_i            = rows;
1213d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1214b1683b59SHong Zhang 
1215d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1216d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1217d2853540SHong Zhang   colorforcol[0] = 0;
1218d2853540SHong Zhang   columns_i      = columns;
1219d2853540SHong Zhang 
1220*e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1221b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1222b1683b59SHong Zhang 
1223b28d80bdSHong Zhang   cm = c->m;
1224b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1225*e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1226b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1227b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1228b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1229b1683b59SHong Zhang     c->ncolumns[i] = n;
1230b1683b59SHong Zhang     if (n) {
1231d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1232b1683b59SHong Zhang     }
1233d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1234d2853540SHong Zhang     columns_i       += n;
1235b1683b59SHong Zhang 
1236b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1237e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1238*e99be685SHong Zhang 
1239b1683b59SHong Zhang     /* loop over columns*/
1240b1683b59SHong Zhang     for (j=0; j<n; j++) {
1241b1683b59SHong Zhang       col     = is[j];
1242d2853540SHong Zhang       row_idx = cj + ci[col];
1243b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1244b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1245b1683b59SHong Zhang       for (k=0; k<m; k++) {
1246*e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1247d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1248b1683b59SHong Zhang       }
1249b1683b59SHong Zhang     }
1250b1683b59SHong Zhang     /* count the number of hits */
1251b1683b59SHong Zhang     nrows = 0;
1252e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1253b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1254b1683b59SHong Zhang     }
1255b1683b59SHong Zhang     c->nrows[i]      = nrows;
1256d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1257d2853540SHong Zhang 
1258b1683b59SHong Zhang     nrows       = 0;
1259e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1260b1683b59SHong Zhang       if (rowhit[j]) {
1261d2853540SHong Zhang         rows_i[nrows]            = j;
1262*e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1263b1683b59SHong Zhang         nrows++;
1264b1683b59SHong Zhang       }
1265b1683b59SHong Zhang     }
1266b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1267d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1268b1683b59SHong Zhang   }
1269*e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1270b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1271b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1272d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1273d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1274d2853540SHong Zhang #endif
1275b28d80bdSHong Zhang 
1276d2853540SHong Zhang   c->colorforrow     = colorforrow;
1277d2853540SHong Zhang   c->rows            = rows;
1278d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1279d2853540SHong Zhang   c->colorforcol     = colorforcol;
1280d2853540SHong Zhang   c->columns         = columns;
1281f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1282b1683b59SHong Zhang   PetscFunctionReturn(0);
1283b1683b59SHong Zhang }
1284