xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision e9e4536c4c589450f5ca315e85fba36ecafcdffb)
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;
205c66b693SKris Buschelman 
215c66b693SKris Buschelman   PetscFunctionBegin;
2226be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
23598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2426be0446SHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
25598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2626be0446SHong Zhang   }
27598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
2801e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
29598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
305c66b693SKris Buschelman   PetscFunctionReturn(0);
315c66b693SKris Buschelman }
32ec01deb9SMatthew Knepley EXTERN_C_END
331c24bd37SHong Zhang 
34e9e4536cSHong Zhang /*
35e9e4536cSHong Zhang  MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ - Get symbolic structure of C=A*B
36e9e4536cSHong Zhang   Input Parameter:
37e9e4536cSHong Zhang .    am, Ai, Aj - number of rows and structure of A
38e9e4536cSHong Zhang .    bm, bn, Bi, Bj - number of rows, columns, and structure of B
39e9e4536cSHong Zhang .    fill - filll ratio See MatMatMult()
40e9e4536cSHong Zhang 
41e9e4536cSHong Zhang   Output Parameter:
42e9e4536cSHong Zhang .    Ci, Cj - structure of C = A*B
43e9e4536cSHong Zhang .    nspacedouble - number of extra mallocs
44e9e4536cSHong Zhang  */
4526be0446SHong Zhang #undef __FUNCT__
46b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ"
47b561aa9dSHong 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)
4858c24d83SHong Zhang {
49dfbe8321SBarry Smith   PetscErrorCode     ierr;
50a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
51b561aa9dSHong Zhang   PetscInt           *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj;
52b561aa9dSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0;
53be0fcf8dSHong Zhang   PetscBT            lnkbt;
5458c24d83SHong Zhang 
5558c24d83SHong Zhang   PetscFunctionBegin;
5658c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
5758c24d83SHong Zhang   /* free space for accumulating nonzero column info */
5838baddfdSBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5958c24d83SHong Zhang   ci[0] = 0;
6058c24d83SHong Zhang 
61be0fcf8dSHong Zhang   /* create and initialize a linked list */
62be0fcf8dSHong Zhang   nlnk = bn+1;
63be0fcf8dSHong Zhang   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
6458c24d83SHong Zhang 
65c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
66a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6758c24d83SHong Zhang   current_space = free_space;
6858c24d83SHong Zhang 
6958c24d83SHong Zhang   /* Determine symbolic info for each row of the product: */
7058c24d83SHong Zhang   for (i=0; i<am; i++) {
7158c24d83SHong Zhang     anzi = ai[i+1] - ai[i];
7258c24d83SHong Zhang     cnzi = 0;
73b561aa9dSHong Zhang     aj   = Aj + ai[i];
7477584fe6SHong Zhang     for (j=0; j<anzi; j++){
75b561aa9dSHong Zhang       brow = aj[j];
7658c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
7758c24d83SHong Zhang       bjj  = bj + bi[brow];
781c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
79dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
801c239cc6SHong Zhang       cnzi += nlnk;
8158c24d83SHong Zhang     }
8258c24d83SHong Zhang 
8358c24d83SHong Zhang     /* If free space is not available, make more free space */
8458c24d83SHong Zhang     /* Double the amount of total space in the list */
8558c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
864238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
87b561aa9dSHong Zhang       ndouble++;
8858c24d83SHong Zhang     }
8958c24d83SHong Zhang 
90c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
91be0fcf8dSHong Zhang     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
92c5db241fSHong Zhang     current_space->array           += cnzi;
9358c24d83SHong Zhang     current_space->local_used      += cnzi;
9458c24d83SHong Zhang     current_space->local_remaining -= cnzi;
9558c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
9658c24d83SHong Zhang   }
9758c24d83SHong Zhang 
9858c24d83SHong Zhang   /* Column indices are in the list of free space */
9958c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
10058c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
10138baddfdSBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
102a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
103be0fcf8dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
10458c24d83SHong Zhang 
105b561aa9dSHong Zhang   *Ci           = ci;
106b561aa9dSHong Zhang   *Cj           = cj;
107b561aa9dSHong Zhang   *nspacedouble = ndouble;
108b561aa9dSHong Zhang   PetscFunctionReturn(0);
109b561aa9dSHong Zhang }
110b561aa9dSHong Zhang 
111b561aa9dSHong Zhang #undef __FUNCT__
112e9e4536cSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_SparseAxpy"
113e9e4536cSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_SparseAxpy(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
114e9e4536cSHong Zhang {
115e9e4536cSHong Zhang   PetscErrorCode ierr;
116e9e4536cSHong Zhang   PetscInt       *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj,rmax=0,*abj,*cj_tmp,nextabj;
117e9e4536cSHong Zhang   PetscInt       i,j,anzi,brow,bnzj,cnzi,k;
118e9e4536cSHong Zhang   PetscBT        bt;
119e9e4536cSHong Zhang 
120e9e4536cSHong Zhang   PetscFunctionBegin;
121e9e4536cSHong Zhang   /* Allocate ci array, arrays for fill computation and */
122e9e4536cSHong Zhang   /* free space for accumulating nonzero column info */
123e9e4536cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
124e9e4536cSHong Zhang   ci[0] = 0;
125e9e4536cSHong Zhang 
126e9e4536cSHong Zhang   /* Get ci and rmax for C */
127e9e4536cSHong Zhang   ierr = PetscBTCreate(bn,bt);CHKERRQ(ierr);
128e9e4536cSHong Zhang   ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
129e9e4536cSHong Zhang   for (i=0; i<am; i++) {
130e9e4536cSHong Zhang     anzi = ai[i+1] - ai[i];
131e9e4536cSHong Zhang     cnzi = 0;
132e9e4536cSHong Zhang     aj   = Aj + ai[i];
133e9e4536cSHong Zhang     for (j=0; j<anzi; j++){
134e9e4536cSHong Zhang       brow = aj[j];
135e9e4536cSHong Zhang       bnzj = bi[brow+1] - bi[brow];
136e9e4536cSHong Zhang       bjj  = bj + bi[brow];
137e9e4536cSHong Zhang       for (k=0; k<bnzj; k++){
138e9e4536cSHong Zhang         if (!PetscBTLookupSet(bt,bjj[k])){  /* new entry */
139e9e4536cSHong Zhang           cnzi++;
140e9e4536cSHong Zhang         }
141e9e4536cSHong Zhang       }
142e9e4536cSHong Zhang     }
143e9e4536cSHong Zhang     ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); /* optimize this? */
144e9e4536cSHong Zhang     ci[i+1] = ci[i] + cnzi;
145e9e4536cSHong Zhang     if (rmax < cnzi) rmax = cnzi;
146e9e4536cSHong Zhang   }
147e9e4536cSHong Zhang 
148e9e4536cSHong Zhang   /* Allocate space for cj */
149e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
150e9e4536cSHong Zhang 
151e9e4536cSHong Zhang   /* allocate a temp array for storing column indices of A*B */
152e9e4536cSHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscInt),&abj);CHKERRQ(ierr);
153e9e4536cSHong Zhang 
154e9e4536cSHong Zhang   /* Determine cj */
155e9e4536cSHong Zhang   for (i=0; i<am; i++) {
156e9e4536cSHong Zhang     anzi = ai[i+1] - ai[i];
157e9e4536cSHong Zhang     cnzi = 0;
158e9e4536cSHong Zhang     nextabj=0;
159e9e4536cSHong Zhang     aj   = Aj + ai[i];
160e9e4536cSHong Zhang     for (j=0; j<anzi; j++){
161e9e4536cSHong Zhang       brow = aj[j];
162e9e4536cSHong Zhang       bnzj = bi[brow+1] - bi[brow];
163e9e4536cSHong Zhang       bjj  = bj + bi[brow];
164e9e4536cSHong Zhang       for (k=0; k<bnzj; k++){
165e9e4536cSHong Zhang         if (!PetscBTLookupSet(bt,bjj[k])){  /* new entry */
166e9e4536cSHong Zhang           abj[nextabj] = bjj[k]; nextabj++;
167e9e4536cSHong Zhang         }
168e9e4536cSHong Zhang       }
169e9e4536cSHong Zhang     }
170e9e4536cSHong Zhang 
171e9e4536cSHong Zhang     /* sort abj, then copy it to cj */
172e9e4536cSHong Zhang     cnzi = ci[i+1] - ci[i];
173e9e4536cSHong Zhang     ierr = PetscSortInt(cnzi,abj);CHKERRQ(ierr);
174e9e4536cSHong Zhang     ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
175e9e4536cSHong Zhang     cj_tmp = cj + ci[i];
176e9e4536cSHong Zhang     for (k=0; k< cnzi; k++){
177e9e4536cSHong Zhang       cj_tmp[k] = abj[k];
178e9e4536cSHong Zhang     }
179e9e4536cSHong Zhang   }
180e9e4536cSHong Zhang 
181e9e4536cSHong Zhang   ierr = PetscBTDestroy(bt);CHKERRQ(ierr);
182e9e4536cSHong Zhang   ierr = PetscFree(abj);CHKERRQ(ierr);
183e9e4536cSHong Zhang 
184e9e4536cSHong Zhang   *Ci           = ci;
185e9e4536cSHong Zhang   *Cj           = cj;
186e9e4536cSHong Zhang   *nspacedouble = 0;
187e9e4536cSHong Zhang   PetscFunctionReturn(0);
188e9e4536cSHong Zhang }
189e9e4536cSHong Zhang 
190e9e4536cSHong Zhang #undef __FUNCT__
191b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
192b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
193b561aa9dSHong Zhang {
194b561aa9dSHong Zhang   PetscErrorCode ierr;
195b561aa9dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
196b561aa9dSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj;
197b561aa9dSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
198b561aa9dSHong Zhang   MatScalar      *ca;
199b561aa9dSHong Zhang   PetscReal      afill;
20036ec6d2dSHong Zhang   PetscBool      dense_axpy; /* false: use sparse axpy; otherwise use dense axpy in MatMatMultNumeric_SeqAIJ_SeqAIJ() */
201b561aa9dSHong Zhang 
202b561aa9dSHong Zhang   PetscFunctionBegin;
203b561aa9dSHong Zhang   /* Get ci and cj */
204b561aa9dSHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
205e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
206e9e4536cSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ() is done \n");CHKERRQ(ierr);
207e9e4536cSHong Zhang #endif
208b561aa9dSHong Zhang 
20958c24d83SHong Zhang   /* Allocate space for ca */
21058c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
21158c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
21258c24d83SHong Zhang 
21326be0446SHong Zhang   /* put together the new symbolic matrix */
2147adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
21558c24d83SHong Zhang 
21658c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
21758c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
21858c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
219e6b907acSBarry Smith   c->free_a   = PETSC_TRUE;
220e6b907acSBarry Smith   c->free_ij  = PETSC_TRUE;
22158c24d83SHong Zhang   c->nonew    = 0;
2228cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ;
22358c24d83SHong Zhang 
22453565b12SHong Zhang   /* Determine which MatMatMultNumeric_SeqAIJ_SeqAIJ() to be used */
22536ec6d2dSHong Zhang   dense_axpy = PETSC_TRUE;
22636ec6d2dSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_denseaxpy",&dense_axpy,PETSC_NULL);CHKERRQ(ierr);
22736ec6d2dSHong Zhang   if (dense_axpy){
22836ec6d2dSHong Zhang     ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr);
229c58ca2e3SHong Zhang     ierr = PetscMemzero(c->matmult_abdense,dense_axpy*bn*sizeof(PetscScalar));CHKERRQ(ierr);
230c58ca2e3SHong Zhang     (*C)->ops->matmultnumeric =  MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */
231c58ca2e3SHong Zhang   } else { /* slower, but use less memory */
232e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
233e9e4536cSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"call  MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy \n");
234e9e4536cSHong Zhang #endif
235c58ca2e3SHong Zhang     (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* slower, less memory */
236c58ca2e3SHong Zhang   }
2370b7e3e3dSHong Zhang 
2387212da7cSHong Zhang   /* set MatInfo */
2397212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
240f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
2417212da7cSHong Zhang   c->maxnz                     = ci[am];
2427212da7cSHong Zhang   c->nz                        = ci[am];
2437212da7cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
2447212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
2457212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
2467212da7cSHong Zhang 
2477212da7cSHong Zhang #if defined(PETSC_USE_INFO)
2487212da7cSHong Zhang   if (ci[am]) {
249f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
250f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
251f2b054eeSHong Zhang   } else {
252f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
253be0fcf8dSHong Zhang   }
254f2b054eeSHong Zhang #endif
25558c24d83SHong Zhang   PetscFunctionReturn(0);
25658c24d83SHong Zhang }
257d50806bdSBarry Smith 
25826be0446SHong Zhang #undef __FUNCT__
25926be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
260dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
261d50806bdSBarry Smith {
262dfbe8321SBarry Smith   PetscErrorCode ierr;
263d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
264d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
265d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
266d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
26738baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
268c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
269519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
27036ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
2710b7e3e3dSHong Zhang   PetscScalar    *ab_dense=c->matmult_abdense;
272d50806bdSBarry Smith 
273d50806bdSBarry Smith   PetscFunctionBegin;
274c124e916SHong Zhang   /* clean old values in C */
275c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
276d50806bdSBarry Smith   /* Traverse A row-wise. */
277d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
278d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
279d50806bdSBarry Smith   for (i=0; i<am; i++) {
280d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
281d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
282519eb980SHong Zhang       brow = aj[j];
283d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
284d50806bdSBarry Smith       bjj  = bj + bi[brow];
285d50806bdSBarry Smith       baj  = ba + bi[brow];
286519eb980SHong Zhang       /* perform dense axpy */
28736ec6d2dSHong Zhang       valtmp = aa[j];
288519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
28936ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
290519eb980SHong Zhang       }
291519eb980SHong Zhang       flops += 2*bnzi;
292519eb980SHong Zhang     }
293c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
294c58ca2e3SHong Zhang 
295c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
296519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
297519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
298519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
299519eb980SHong Zhang     }
300519eb980SHong Zhang     flops += cnzi;
301519eb980SHong Zhang     cj += cnzi; ca += cnzi;
302519eb980SHong Zhang   }
303c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
304c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
305c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
306c58ca2e3SHong Zhang   C->ops->matmultnumeric =  MatMatMultNumeric_SeqAIJ_SeqAIJ;
307c58ca2e3SHong Zhang   PetscFunctionReturn(0);
308c58ca2e3SHong Zhang }
309c58ca2e3SHong Zhang 
310c58ca2e3SHong Zhang #undef __FUNCT__
311c58ca2e3SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
312c58ca2e3SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat B,Mat C)
313c58ca2e3SHong Zhang {
314c58ca2e3SHong Zhang   PetscErrorCode ierr;
315c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
316c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
317c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
318c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
319c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
320c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
321c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
32236ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
323c58ca2e3SHong Zhang   PetscInt       nextb;
324c58ca2e3SHong Zhang 
325c58ca2e3SHong Zhang   PetscFunctionBegin;
326e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
327e9e4536cSHong Zhang   //ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy...\n");CHKERRQ(ierr);
328e9e4536cSHong Zhang #endif
329c58ca2e3SHong Zhang   /* clean old values in C */
330c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
331c58ca2e3SHong Zhang   /* Traverse A row-wise. */
332c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
333c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
334519eb980SHong Zhang   for (i=0;i<am;i++) {
335519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
336519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
337519eb980SHong Zhang     for (j=0;j<anzi;j++) {
338519eb980SHong Zhang       brow = aj[j];
339519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
340519eb980SHong Zhang       bjj  = bj + bi[brow];
341519eb980SHong Zhang       baj  = ba + bi[brow];
342519eb980SHong Zhang       /* perform sparse axpy */
34336ec6d2dSHong Zhang       valtmp = aa[j];
344c124e916SHong Zhang       nextb  = 0;
345c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
346c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
34736ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
348c124e916SHong Zhang         }
349d50806bdSBarry Smith       }
350d50806bdSBarry Smith       flops += 2*bnzi;
351d50806bdSBarry Smith     }
352519eb980SHong Zhang     aj += anzi; aa += anzi;
353519eb980SHong Zhang     cj += cnzi; ca += cnzi;
354d50806bdSBarry Smith   }
355c58ca2e3SHong Zhang 
356716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
357716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
358d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
359d50806bdSBarry Smith   PetscFunctionReturn(0);
360d50806bdSBarry Smith }
361bc011b1eSHong Zhang 
362e9e4536cSHong Zhang #undef __FUNCT__
363e9e4536cSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_SparseAxpy"
364e9e4536cSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat B,PetscReal fill,Mat *C)
365e9e4536cSHong Zhang {
366e9e4536cSHong Zhang   PetscErrorCode ierr;
367e9e4536cSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
368e9e4536cSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj;
369e9e4536cSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
370e9e4536cSHong Zhang   MatScalar      *ca;
371e9e4536cSHong Zhang   PetscReal      afill;
372e9e4536cSHong Zhang 
373e9e4536cSHong Zhang   PetscFunctionBegin;
374e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
375e9e4536cSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultSymbolic_SeqAIJ_SeqAIJ_SparseAxpy \n");CHKERRQ(ierr);
376e9e4536cSHong Zhang #endif
377e9e4536cSHong Zhang   /* Get ci and cj */
378e9e4536cSHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_SparseAxpy(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
379e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
380e9e4536cSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_SparseAxpy() is done \n");CHKERRQ(ierr);
381e9e4536cSHong Zhang #endif
382e9e4536cSHong Zhang 
383e9e4536cSHong Zhang   /* Allocate space for ca */
384e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
385e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
386e9e4536cSHong Zhang 
387e9e4536cSHong Zhang   /* put together the new symbolic matrix */
388e9e4536cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
389e9e4536cSHong Zhang 
390e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
391e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
392e9e4536cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
393e9e4536cSHong Zhang   c->free_a   = PETSC_TRUE;
394e9e4536cSHong Zhang   c->free_ij  = PETSC_TRUE;
395e9e4536cSHong Zhang   c->nonew    = 0;
396e9e4536cSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* slower, less memory */
397e9e4536cSHong Zhang 
398e9e4536cSHong Zhang   /* set MatInfo */
399e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
400e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
401e9e4536cSHong Zhang   c->maxnz                     = ci[am];
402e9e4536cSHong Zhang   c->nz                        = ci[am];
403e9e4536cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
404e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
405e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
406e9e4536cSHong Zhang 
407e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
408e9e4536cSHong Zhang   if (ci[am]) {
409e9e4536cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
410e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
411e9e4536cSHong Zhang   } else {
412e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
413e9e4536cSHong Zhang   }
414e9e4536cSHong Zhang #endif
415e9e4536cSHong Zhang   PetscFunctionReturn(0);
416e9e4536cSHong Zhang }
417e9e4536cSHong Zhang 
418e9e4536cSHong Zhang 
419d2853540SHong Zhang /* This routine is not used. Should be removed! */
420bc011b1eSHong Zhang #undef __FUNCT__
4216fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
4226fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4235df89d91SHong Zhang {
424bc011b1eSHong Zhang   PetscErrorCode ierr;
425bc011b1eSHong Zhang 
426bc011b1eSHong Zhang   PetscFunctionBegin;
427bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
4286fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
429bc011b1eSHong Zhang   }
4306fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
431bc011b1eSHong Zhang   PetscFunctionReturn(0);
432bc011b1eSHong Zhang }
433bc011b1eSHong Zhang 
434bc011b1eSHong Zhang #undef __FUNCT__
435e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
436e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
4372128a86cSHong Zhang {
4382128a86cSHong Zhang   PetscErrorCode      ierr;
439e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
4402128a86cSHong Zhang 
4412128a86cSHong Zhang   PetscFunctionBegin;
442b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
4432128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
4442128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
4452128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
4462128a86cSHong Zhang   PetscFunctionReturn(0);
4472128a86cSHong Zhang }
4482128a86cSHong Zhang 
4492128a86cSHong Zhang #undef __FUNCT__
4502128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
4512128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
4522128a86cSHong Zhang {
4532128a86cSHong Zhang   PetscErrorCode      ierr;
4542128a86cSHong Zhang   PetscContainer      container;
455e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
4562128a86cSHong Zhang 
4572128a86cSHong Zhang   PetscFunctionBegin;
458e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
4592128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
4602128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
4612128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
4622128a86cSHong Zhang   if (A->ops->destroy) {
4632128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
4642128a86cSHong Zhang   }
465e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
4662128a86cSHong Zhang   PetscFunctionReturn(0);
4672128a86cSHong Zhang }
4682128a86cSHong Zhang 
4692128a86cSHong Zhang #undef __FUNCT__
4706fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
4716fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
472bc011b1eSHong Zhang {
4735df89d91SHong Zhang   PetscErrorCode      ierr;
47481d82fe4SHong Zhang   Mat                 Bt;
47581d82fe4SHong Zhang   PetscInt            *bti,*btj;
476e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
4772128a86cSHong Zhang   PetscContainer      container;
478d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
479d2853540SHong Zhang 
48081d82fe4SHong Zhang   PetscFunctionBegin;
481d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
48281d82fe4SHong Zhang    /* create symbolic Bt */
48381d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
48481d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
48581d82fe4SHong Zhang 
48681d82fe4SHong Zhang   /* get symbolic C=A*Bt */
48781d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
48881d82fe4SHong Zhang 
4892128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
490e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
4912128a86cSHong Zhang 
4922128a86cSHong Zhang   /* attach the supporting struct to C */
4932128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
4942128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
495e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
496e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
4972128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4982128a86cSHong Zhang 
4992128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
5002128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
5012128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
5022128a86cSHong Zhang 
503d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
504d2853540SHong Zhang   etime2 += tf - t0;
505d2853540SHong Zhang 
506f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
5072128a86cSHong Zhang   if (multtrans->usecoloring){
508b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
509b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
5102128a86cSHong Zhang     ISColoring           iscoloring;
5112128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
512d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
513e8354b3eSHong Zhang 
514e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
515502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
516d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
517d2853540SHong Zhang     etime0 += tf - t0;
518d2853540SHong Zhang 
519d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
520b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
5212128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
5222128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
523e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
524d2853540SHong Zhang     etime01 += tf - t0;
5252128a86cSHong Zhang 
526e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
5272128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
5282128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
5292128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
5302128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
5312128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
5322128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
5332128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
5342128a86cSHong Zhang 
5352128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
5362128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
5372128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
5382128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
5392128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
5402128a86cSHong Zhang     multtrans->ABt_den = C_dense;
541e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
542e8354b3eSHong Zhang     etime1 += tf - t0;
543f94ccd6cSHong Zhang 
544f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
5451ce71dffSSatish Balay     {
546f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
547f94ccd6cSHong 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));
548f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
5491ce71dffSSatish Balay     }
550f94ccd6cSHong Zhang #endif
5512128a86cSHong Zhang   }
552*e99be685SHong Zhang   /* clean up */
553*e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
554*e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
5552128a86cSHong Zhang 
556f94ccd6cSHong Zhang 
557f94ccd6cSHong Zhang 
55881d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
55981d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
5601fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
5611fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5621fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
5631fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
5641fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
5651fa1209cSHong Zhang   MatScalar          *ca;
5661fa1209cSHong Zhang   PetscBT            lnkbt;
5671fa1209cSHong Zhang   PetscReal          afill;
5685df89d91SHong Zhang 
5691fa1209cSHong Zhang   /* Allocate row pointer array ci  */
5701fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5711fa1209cSHong Zhang   ci[0] = 0;
5721fa1209cSHong Zhang 
5731fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
5741fa1209cSHong Zhang   nlnk = bm+1;
5751fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5761fa1209cSHong Zhang 
5771fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
5781fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5791fa1209cSHong Zhang   current_space = free_space;
5801fa1209cSHong Zhang 
5811fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
5821fa1209cSHong Zhang   for (i=0; i<am; i++) {
5831fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
5841fa1209cSHong Zhang     cnzi = 0;
5851fa1209cSHong Zhang     acol = aj + ai[i];
5861fa1209cSHong Zhang     for (j=0; j<bm; j++){
5871fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
5881fa1209cSHong Zhang       bcol= bj + bi[j];
58981d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
5901fa1209cSHong Zhang       ka = 0; kb = 0;
5911fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
59281d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
59381d82fe4SHong Zhang         if (ka == anzi) break;
59481d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
59581d82fe4SHong Zhang         if (kb == bnzj) break;
5961fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
5971fa1209cSHong Zhang           index[0] = j;
5981fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5991fa1209cSHong Zhang           cnzi++;
6001fa1209cSHong Zhang           break;
6011fa1209cSHong Zhang         }
6021fa1209cSHong Zhang       }
6031fa1209cSHong Zhang     }
6041fa1209cSHong Zhang 
6051fa1209cSHong Zhang     /* If free space is not available, make more free space */
6061fa1209cSHong Zhang     /* Double the amount of total space in the list */
6071fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
6081fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
6091fa1209cSHong Zhang       nspacedouble++;
6101fa1209cSHong Zhang     }
6111fa1209cSHong Zhang 
6121fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
6131fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
6141fa1209cSHong Zhang     current_space->array           += cnzi;
6151fa1209cSHong Zhang     current_space->local_used      += cnzi;
6161fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
6171fa1209cSHong Zhang 
6181fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
6191fa1209cSHong Zhang   }
6201fa1209cSHong Zhang 
6211fa1209cSHong Zhang 
6221fa1209cSHong Zhang   /* Column indices are in the list of free space.
6231fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
6241fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6251fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6261fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
6271fa1209cSHong Zhang 
6281fa1209cSHong Zhang   /* Allocate space for ca */
6291fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
6301fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
6311fa1209cSHong Zhang 
6321fa1209cSHong Zhang   /* put together the new symbolic matrix */
6331fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
6341fa1209cSHong Zhang 
6351fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6361fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6371fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
6381fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
6391fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
6401fa1209cSHong Zhang   c->nonew    = 0;
6411fa1209cSHong Zhang 
6421fa1209cSHong Zhang   /* set MatInfo */
6431fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6441fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
6451fa1209cSHong Zhang   c->maxnz                     = ci[am];
6461fa1209cSHong Zhang   c->nz                        = ci[am];
6471fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
6481fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
6491fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
6501fa1209cSHong Zhang 
6511fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
6521fa1209cSHong Zhang   if (ci[am]) {
6531fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
6546fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
6551fa1209cSHong Zhang   } else {
6561fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6571fa1209cSHong Zhang   }
6581fa1209cSHong Zhang #endif
65981d82fe4SHong Zhang #endif
6605df89d91SHong Zhang   PetscFunctionReturn(0);
6615df89d91SHong Zhang }
6625df89d91SHong Zhang 
6636973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
6645df89d91SHong Zhang #undef __FUNCT__
6656fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
6666fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
6675df89d91SHong Zhang {
6685df89d91SHong Zhang   PetscErrorCode ierr;
6695df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
670e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
67181d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
6725df89d91SHong Zhang   PetscLogDouble flops=0.0;
6736973769bSHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
674e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
6752128a86cSHong Zhang   PetscContainer      container;
6766973769bSHong Zhang #if defined(USE_ARRAY)
6776973769bSHong Zhang   MatScalar      *spdot;
6786973769bSHong Zhang #endif
6795df89d91SHong Zhang 
6805df89d91SHong Zhang   PetscFunctionBegin;
681e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
6822128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
6832128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
6842128a86cSHong Zhang   if (multtrans->usecoloring){
685b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
686c8db22f6SHong Zhang     Mat                   Bt_dense;
687c8db22f6SHong Zhang     PetscInt              m,n;
688b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
689a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
690a0b95ffbSSatish Balay 
6912128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
692c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
693c8db22f6SHong Zhang 
694b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
6952128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
696b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
6972128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
698b2d2b04fSHong Zhang     etime0 += tf - t0;
699c8db22f6SHong Zhang 
700c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
7012128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
7022128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
7032128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7042128a86cSHong Zhang     etime2 += tf - t0;
705c8db22f6SHong Zhang 
7062128a86cSHong Zhang     /* Recover C from C_dense */
7072128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
708b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
7092128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7102128a86cSHong Zhang     etime1 += tf - t0;
711f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
712f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
713f94ccd6cSHong Zhang #endif
714c8db22f6SHong Zhang     PetscFunctionReturn(0);
715c8db22f6SHong Zhang   }
716c8db22f6SHong Zhang 
7176973769bSHong Zhang #if defined(USE_ARRAY)
7186973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
719e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
720e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
7216973769bSHong Zhang #endif
7226973769bSHong Zhang 
72381d82fe4SHong Zhang   /* clear old values in C */
72481d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
7255df89d91SHong Zhang 
7261fa1209cSHong Zhang   for (i=0; i<cm; i++) {
72781d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
72881d82fe4SHong Zhang     acol = aj + ai[i];
7296973769bSHong Zhang     aval = aa + ai[i];
7301fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
7311fa1209cSHong Zhang     ccol = cj + ci[i];
7326973769bSHong Zhang     cval = ca + ci[i];
7331fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
73481d82fe4SHong Zhang       brow = ccol[j];
73581d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
73681d82fe4SHong Zhang       bcol = bj + bi[brow];
7376973769bSHong Zhang       bval = ba + bi[brow];
7386973769bSHong Zhang 
73981d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
7406973769bSHong Zhang #if defined(USE_ARRAY)
7416973769bSHong Zhang       /* put ba to spdot array */
7426973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
7436973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
7446973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
7456973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
7466973769bSHong Zhang       }
7476973769bSHong Zhang       /* zero spdot array */
7486973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
7496973769bSHong Zhang #else
75081d82fe4SHong Zhang       nexta = 0; nextb = 0;
75181d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
75281d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
75381d82fe4SHong Zhang         if (nexta == anzi) break;
75481d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
75581d82fe4SHong Zhang         if (nextb == bnzj) break;
75681d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
7576973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
75881d82fe4SHong Zhang           nexta++; nextb++;
75981d82fe4SHong Zhang           flops += 2;
7601fa1209cSHong Zhang         }
7611fa1209cSHong Zhang       }
7626973769bSHong Zhang #endif
76381d82fe4SHong Zhang     }
76481d82fe4SHong Zhang   }
76581d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76681d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76781d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
7686973769bSHong Zhang #if defined(USE_ARRAY)
7696973769bSHong Zhang   ierr = PetscFree(spdot);
7706973769bSHong Zhang #endif
7715df89d91SHong Zhang   PetscFunctionReturn(0);
7725df89d91SHong Zhang }
7735df89d91SHong Zhang 
7745df89d91SHong Zhang #undef __FUNCT__
77575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
77675648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
7775df89d91SHong Zhang   PetscErrorCode ierr;
7785df89d91SHong Zhang 
7795df89d91SHong Zhang   PetscFunctionBegin;
7805df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
78175648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
7825df89d91SHong Zhang   }
78375648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
7845df89d91SHong Zhang   PetscFunctionReturn(0);
7855df89d91SHong Zhang }
7865df89d91SHong Zhang 
7875df89d91SHong Zhang #undef __FUNCT__
78875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
78975648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
7905df89d91SHong Zhang {
791bc011b1eSHong Zhang   PetscErrorCode ierr;
792bc011b1eSHong Zhang   Mat            At;
79338baddfdSBarry Smith   PetscInt       *ati,*atj;
794bc011b1eSHong Zhang 
795bc011b1eSHong Zhang   PetscFunctionBegin;
796bc011b1eSHong Zhang   /* create symbolic At */
797bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
798d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
799bc011b1eSHong Zhang 
800bc011b1eSHong Zhang   /* get symbolic C=At*B */
801bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
802bc011b1eSHong Zhang 
803bc011b1eSHong Zhang   /* clean up */
8046bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
805bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
806bc011b1eSHong Zhang   PetscFunctionReturn(0);
807bc011b1eSHong Zhang }
808bc011b1eSHong Zhang 
809bc011b1eSHong Zhang #undef __FUNCT__
81075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
81175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
812bc011b1eSHong Zhang {
813bc011b1eSHong Zhang   PetscErrorCode ierr;
8140fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
815d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
816d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
817d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
8180fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
819bc011b1eSHong Zhang 
820bc011b1eSHong Zhang   PetscFunctionBegin;
821bc011b1eSHong Zhang   /* clear old values in C */
822bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
823bc011b1eSHong Zhang 
824bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
825bc011b1eSHong Zhang   for (i=0;i<am;i++) {
826bc011b1eSHong Zhang     bj   = b->j + bi[i];
827bc011b1eSHong Zhang     ba   = b->a + bi[i];
828bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
829bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
830bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
831bc011b1eSHong Zhang       nextb = 0;
8320fbc74f4SHong Zhang       crow  = *aj++;
833bc011b1eSHong Zhang       cjj   = cj + ci[crow];
834bc011b1eSHong Zhang       caj   = ca + ci[crow];
835bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
836bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
8370fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
8380fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
839bc011b1eSHong Zhang           nextb++;
840bc011b1eSHong Zhang         }
841bc011b1eSHong Zhang       }
842bc011b1eSHong Zhang       flops += 2*bnzi;
8430fbc74f4SHong Zhang       aa++;
844bc011b1eSHong Zhang     }
845bc011b1eSHong Zhang   }
846bc011b1eSHong Zhang 
847bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
848bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
849bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
850bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
851bc011b1eSHong Zhang   PetscFunctionReturn(0);
852bc011b1eSHong Zhang }
8539123193aSHong Zhang 
854ec01deb9SMatthew Knepley EXTERN_C_BEGIN
8559123193aSHong Zhang #undef __FUNCT__
8569123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
8579123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8589123193aSHong Zhang {
8599123193aSHong Zhang   PetscErrorCode ierr;
8609123193aSHong Zhang 
8619123193aSHong Zhang   PetscFunctionBegin;
8629123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
8639123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
8649123193aSHong Zhang   }
8659123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
8669123193aSHong Zhang   PetscFunctionReturn(0);
8679123193aSHong Zhang }
868ec01deb9SMatthew Knepley EXTERN_C_END
869edd81eecSMatthew Knepley 
8709123193aSHong Zhang #undef __FUNCT__
8719123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
8729123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
8739123193aSHong Zhang {
8749123193aSHong Zhang   PetscErrorCode ierr;
8759123193aSHong Zhang 
8769123193aSHong Zhang   PetscFunctionBegin;
8775a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
8788cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
8799123193aSHong Zhang   PetscFunctionReturn(0);
8809123193aSHong Zhang }
8819123193aSHong Zhang 
8829123193aSHong Zhang #undef __FUNCT__
8839123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
8849123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
8859123193aSHong Zhang {
886f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
8879123193aSHong Zhang   PetscErrorCode ierr;
888dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
889dd6ea824SBarry Smith   MatScalar      *aa;
890d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
891f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
8929123193aSHong Zhang 
8939123193aSHong Zhang   PetscFunctionBegin;
894f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
895e32f2f54SBarry 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);
896e32f2f54SBarry 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);
897e32f2f54SBarry 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);
898f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
899f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
900f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
901f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
902f32d5d43SBarry Smith     colam = col*am;
903f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
904f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
905f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
906f32d5d43SBarry Smith       aj  = a->j + a->i[i];
907f32d5d43SBarry Smith       aa  = a->a + a->i[i];
908f32d5d43SBarry Smith       for (j=0; j<n; j++) {
909f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
910f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
911f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
912f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
9139123193aSHong Zhang       }
914f32d5d43SBarry Smith       c[colam + i]       = r1;
915f32d5d43SBarry Smith       c[colam + am + i]  = r2;
916f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
917f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
918f32d5d43SBarry Smith     }
919f32d5d43SBarry Smith     b1 += bm4;
920f32d5d43SBarry Smith     b2 += bm4;
921f32d5d43SBarry Smith     b3 += bm4;
922f32d5d43SBarry Smith     b4 += bm4;
923f32d5d43SBarry Smith   }
924f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
925f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
926f32d5d43SBarry Smith       r1 = 0.0;
927f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
928f32d5d43SBarry Smith       aj  = a->j + a->i[i];
929f32d5d43SBarry Smith       aa  = a->a + a->i[i];
930f32d5d43SBarry Smith 
931f32d5d43SBarry Smith       for (j=0; j<n; j++) {
932f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
933f32d5d43SBarry Smith       }
934f32d5d43SBarry Smith       c[col*am + i]     = r1;
935f32d5d43SBarry Smith     }
936f32d5d43SBarry Smith     b1 += bm;
937f32d5d43SBarry Smith   }
938dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
939f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
940f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
941f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
942f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
943f32d5d43SBarry Smith   PetscFunctionReturn(0);
944f32d5d43SBarry Smith }
945f32d5d43SBarry Smith 
946f32d5d43SBarry Smith /*
9474324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
948f32d5d43SBarry Smith */
949f32d5d43SBarry Smith #undef __FUNCT__
950f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
951f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
952f32d5d43SBarry Smith {
953f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
954f32d5d43SBarry Smith   PetscErrorCode ierr;
955dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
956dd6ea824SBarry Smith   MatScalar      *aa;
957d0f46423SBarry 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;
9584324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
959f32d5d43SBarry Smith 
960f32d5d43SBarry Smith   PetscFunctionBegin;
961f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
962f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
963f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
964f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
9654324174eSBarry Smith 
9664324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
9674324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
9684324174eSBarry Smith       colam = col*am;
9694324174eSBarry Smith       arm   = a->compressedrow.nrows;
9704324174eSBarry Smith       ii    = a->compressedrow.i;
9714324174eSBarry Smith       ridx  = a->compressedrow.rindex;
9724324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
9734324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
9744324174eSBarry Smith 	n   = ii[i+1] - ii[i];
9754324174eSBarry Smith 	aj  = a->j + ii[i];
9764324174eSBarry Smith 	aa  = a->a + ii[i];
9774324174eSBarry Smith 	for (j=0; j<n; j++) {
9784324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
9794324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
9804324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
9814324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
9824324174eSBarry Smith 	}
9834324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
9844324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
9854324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
9864324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
9874324174eSBarry Smith       }
9884324174eSBarry Smith       b1 += bm4;
9894324174eSBarry Smith       b2 += bm4;
9904324174eSBarry Smith       b3 += bm4;
9914324174eSBarry Smith       b4 += bm4;
9924324174eSBarry Smith     }
9934324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
9944324174eSBarry Smith       colam = col*am;
9954324174eSBarry Smith       arm   = a->compressedrow.nrows;
9964324174eSBarry Smith       ii    = a->compressedrow.i;
9974324174eSBarry Smith       ridx  = a->compressedrow.rindex;
9984324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
9994324174eSBarry Smith 	r1 = 0.0;
10004324174eSBarry Smith 	n   = ii[i+1] - ii[i];
10014324174eSBarry Smith 	aj  = a->j + ii[i];
10024324174eSBarry Smith 	aa  = a->a + ii[i];
10034324174eSBarry Smith 
10044324174eSBarry Smith 	for (j=0; j<n; j++) {
10054324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
10064324174eSBarry Smith 	}
10074324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
10084324174eSBarry Smith       }
10094324174eSBarry Smith       b1 += bm;
10104324174eSBarry Smith     }
10114324174eSBarry Smith   } else {
1012f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
1013f32d5d43SBarry Smith       colam = col*am;
1014f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1015f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
1016f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1017f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1018f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1019f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1020f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
1021f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
1022f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
1023f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
1024f32d5d43SBarry Smith 	}
1025f32d5d43SBarry Smith 	c[colam + i]       += r1;
1026f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
1027f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
1028f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
1029f32d5d43SBarry Smith       }
1030f32d5d43SBarry Smith       b1 += bm4;
1031f32d5d43SBarry Smith       b2 += bm4;
1032f32d5d43SBarry Smith       b3 += bm4;
1033f32d5d43SBarry Smith       b4 += bm4;
1034f32d5d43SBarry Smith     }
1035f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
1036f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1037f32d5d43SBarry Smith 	r1 = 0.0;
1038f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1039f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1040f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1041f32d5d43SBarry Smith 
1042f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1043f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
1044f32d5d43SBarry Smith 	}
1045f32d5d43SBarry Smith 	c[col*am + i]     += r1;
1046f32d5d43SBarry Smith       }
1047f32d5d43SBarry Smith       b1 += bm;
1048f32d5d43SBarry Smith     }
10494324174eSBarry Smith   }
1050dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
1051f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1052f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
10539123193aSHong Zhang   PetscFunctionReturn(0);
10549123193aSHong Zhang }
1055b1683b59SHong Zhang 
1056b1683b59SHong Zhang #undef __FUNCT__
1057b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1058b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1059c8db22f6SHong Zhang {
1060c8db22f6SHong Zhang   PetscErrorCode ierr;
10612f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
10622f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
10632f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
10642f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
10652f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
10662f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1067c8db22f6SHong Zhang 
1068c8db22f6SHong Zhang   PetscFunctionBegin;
10692f5f1f90SHong Zhang   btval_den=btdense->v;
10702f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
10712f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
10722f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
10732f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1074d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
10752f5f1f90SHong Zhang       btcol = bj + bi[col];
10762f5f1f90SHong Zhang       btval = ba + bi[col];
10772f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1078d2853540SHong Zhang       for (j=0; j<anz; j++){
10792f5f1f90SHong Zhang         brow            = btcol[j];
10802f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1081c8db22f6SHong Zhang       }
1082c8db22f6SHong Zhang     }
10832f5f1f90SHong Zhang     btval_den += m;
1084c8db22f6SHong Zhang   }
1085c8db22f6SHong Zhang   PetscFunctionReturn(0);
1086c8db22f6SHong Zhang }
1087c8db22f6SHong Zhang 
1088c8db22f6SHong Zhang #undef __FUNCT__
1089b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1090b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
10918972f759SHong Zhang {
10928972f759SHong Zhang   PetscErrorCode ierr;
1093b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
10942f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1095b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
10962f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
10978972f759SHong Zhang 
10988972f759SHong Zhang   PetscFunctionBegin;
10998972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
1100b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
1101b2d2b04fSHong Zhang   cp_den = ca_den;
11022f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
11032f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
11042f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
11052f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
11062f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
11072f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1108b2d2b04fSHong Zhang     }
1109b2d2b04fSHong Zhang     cp_den += m;
1110b2d2b04fSHong Zhang   }
1111b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
11128972f759SHong Zhang   PetscFunctionReturn(0);
11138972f759SHong Zhang }
11148972f759SHong Zhang 
1115*e99be685SHong Zhang /*
1116*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1117*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1118*e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1119*e99be685SHong Zhang  */
1120*e99be685SHong Zhang #undef __FUNCT__
1121*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
1122*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1123*e99be685SHong Zhang {
1124*e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1125*e99be685SHong Zhang   PetscErrorCode ierr;
1126*e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1127*e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1128*e99be685SHong Zhang   PetscInt       *cspidx;
1129*e99be685SHong Zhang 
1130*e99be685SHong Zhang   PetscFunctionBegin;
1131*e99be685SHong Zhang   *nn = n;
1132*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1133*e99be685SHong Zhang   if (symmetric) {
1134*e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
1135*e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
1136*e99be685SHong Zhang   } else {
1137*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1138*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1139*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1140*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1141*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1142*e99be685SHong Zhang     jj = a->j;
1143*e99be685SHong Zhang     for (i=0; i<nz; i++) {
1144*e99be685SHong Zhang       collengths[jj[i]]++;
1145*e99be685SHong Zhang     }
1146*e99be685SHong Zhang     cia[0] = oshift;
1147*e99be685SHong Zhang     for (i=0; i<n; i++) {
1148*e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1149*e99be685SHong Zhang     }
1150*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1151*e99be685SHong Zhang     jj   = a->j;
1152*e99be685SHong Zhang     for (row=0; row<m; row++) {
1153*e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1154*e99be685SHong Zhang       for (i=0; i<mr; i++) {
1155*e99be685SHong Zhang         col = *jj++;
1156*e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1157*e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1158*e99be685SHong Zhang       }
1159*e99be685SHong Zhang     }
1160*e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1161*e99be685SHong Zhang     *ia = cia; *ja = cja;
1162*e99be685SHong Zhang     *spidx = cspidx;
1163*e99be685SHong Zhang   }
1164*e99be685SHong Zhang   PetscFunctionReturn(0);
1165*e99be685SHong Zhang }
1166*e99be685SHong Zhang 
1167*e99be685SHong Zhang #undef __FUNCT__
1168*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1169*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1170*e99be685SHong Zhang {
1171*e99be685SHong Zhang   PetscErrorCode ierr;
1172*e99be685SHong Zhang 
1173*e99be685SHong Zhang   PetscFunctionBegin;
1174*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1175*e99be685SHong Zhang 
1176*e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1177*e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1178*e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1179*e99be685SHong Zhang   PetscFunctionReturn(0);
1180*e99be685SHong Zhang }
1181*e99be685SHong Zhang 
11828972f759SHong Zhang #undef __FUNCT__
1183b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1184b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1185b1683b59SHong Zhang {
1186b1683b59SHong Zhang   PetscErrorCode ierr;
1187d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1188b1683b59SHong Zhang   const PetscInt *is;
1189b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1190b1683b59SHong Zhang   IS             *isa;
1191b1683b59SHong Zhang   PetscBool      done;
1192b1683b59SHong Zhang   PetscBool      flg1,flg2;
1193b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1194*e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1195d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1196b1683b59SHong Zhang 
1197b1683b59SHong Zhang   PetscFunctionBegin;
1198b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1199*e99be685SHong Zhang 
1200b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1201b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1202b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1203b1683b59SHong Zhang   if (flg1 || flg2) {
1204b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1205b1683b59SHong Zhang   }
1206b1683b59SHong Zhang 
1207b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1208b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1209b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1210b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1211b1683b59SHong Zhang   c->rstart = 0;
1212b1683b59SHong Zhang 
1213b1683b59SHong Zhang   c->ncolors = nis;
1214b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1215b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1216d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1217d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1218d2853540SHong Zhang   colorforrow[0]    = 0;
1219d2853540SHong Zhang   rows_i            = rows;
1220d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1221b1683b59SHong Zhang 
1222d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1223d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1224d2853540SHong Zhang   colorforcol[0] = 0;
1225d2853540SHong Zhang   columns_i      = columns;
1226d2853540SHong Zhang 
1227*e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1228b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1229b1683b59SHong Zhang 
1230b28d80bdSHong Zhang   cm = c->m;
1231b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1232*e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1233b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1234b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1235b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1236b1683b59SHong Zhang     c->ncolumns[i] = n;
1237b1683b59SHong Zhang     if (n) {
1238d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1239b1683b59SHong Zhang     }
1240d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1241d2853540SHong Zhang     columns_i       += n;
1242b1683b59SHong Zhang 
1243b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1244e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1245*e99be685SHong Zhang 
1246b1683b59SHong Zhang     /* loop over columns*/
1247b1683b59SHong Zhang     for (j=0; j<n; j++) {
1248b1683b59SHong Zhang       col     = is[j];
1249d2853540SHong Zhang       row_idx = cj + ci[col];
1250b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1251b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1252b1683b59SHong Zhang       for (k=0; k<m; k++) {
1253*e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1254d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1255b1683b59SHong Zhang       }
1256b1683b59SHong Zhang     }
1257b1683b59SHong Zhang     /* count the number of hits */
1258b1683b59SHong Zhang     nrows = 0;
1259e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1260b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1261b1683b59SHong Zhang     }
1262b1683b59SHong Zhang     c->nrows[i]      = nrows;
1263d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1264d2853540SHong Zhang 
1265b1683b59SHong Zhang     nrows       = 0;
1266e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1267b1683b59SHong Zhang       if (rowhit[j]) {
1268d2853540SHong Zhang         rows_i[nrows]            = j;
1269*e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1270b1683b59SHong Zhang         nrows++;
1271b1683b59SHong Zhang       }
1272b1683b59SHong Zhang     }
1273b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1274d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1275b1683b59SHong Zhang   }
1276*e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1277b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1278b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1279d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1280d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1281d2853540SHong Zhang #endif
1282b28d80bdSHong Zhang 
1283d2853540SHong Zhang   c->colorforrow     = colorforrow;
1284d2853540SHong Zhang   c->rows            = rows;
1285d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1286d2853540SHong Zhang   c->colorforcol     = colorforcol;
1287d2853540SHong Zhang   c->columns         = columns;
1288f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1289b1683b59SHong Zhang   PetscFunctionReturn(0);
1290b1683b59SHong Zhang }
1291