xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 251f4c6745e69799cb0300fd80598c58f96eb452)
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;
203c50cad2SHong Zhang   PetscBool      scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE;
215c66b693SKris Buschelman 
225c66b693SKris Buschelman   PetscFunctionBegin;
2326be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
24152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
25152983bfSHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
263c50cad2SHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,PETSC_NULL);CHKERRQ(ierr);
27d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
28d8bbc50fSBarry Smith     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
293c50cad2SHong Zhang     if (scalable_fast){
303c50cad2SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
31aacf854cSBarry Smith     } else if (scalable){
32aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
3325023028SHong Zhang     } else {
3426be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3525023028SHong Zhang     }
36598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3726be0446SHong Zhang   }
38598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3901e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
40598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
415c66b693SKris Buschelman   PetscFunctionReturn(0);
425c66b693SKris Buschelman }
43ec01deb9SMatthew Knepley EXTERN_C_END
441c24bd37SHong Zhang 
45e9e4536cSHong Zhang #undef __FUNCT__
46b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
47b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
48b561aa9dSHong Zhang {
49b561aa9dSHong Zhang   PetscErrorCode     ierr;
50b561aa9dSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
51971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
5225c41797SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
53b561aa9dSHong Zhang   PetscReal          afill;
54fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
55fb908031SHong Zhang   PetscBT            lnkbt;
56fb908031SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
57b561aa9dSHong Zhang 
58b561aa9dSHong Zhang   PetscFunctionBegin;
59bd958071SHong Zhang   /* Get ci and cj */
60bd958071SHong Zhang   /*---------------*/
61fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
62fb908031SHong Zhang   /* free space for accumulating nonzero column info */
63fb908031SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
64fb908031SHong Zhang   ci[0] = 0;
65fb908031SHong Zhang 
66fb908031SHong Zhang   /* create and initialize a linked list */
67fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
68fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
69fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
70fb908031SHong Zhang 
71fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
72fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
73fb908031SHong Zhang   current_space = free_space;
74fb908031SHong Zhang 
75fb908031SHong Zhang   /* Determine ci and cj */
76fb908031SHong Zhang   for (i=0; i<am; i++) {
77fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
78fb908031SHong Zhang     aj   = a->j + ai[i];
79fb908031SHong Zhang     for (j=0; j<anzi; j++){
80fb908031SHong Zhang       brow = aj[j];
81fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
82fb908031SHong Zhang       bj   = b->j + bi[brow];
83fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
84fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
85fb908031SHong Zhang     }
86fb908031SHong Zhang     cnzi = lnk[0];
87fb908031SHong Zhang 
88fb908031SHong Zhang     /* If free space is not available, make more free space */
89fb908031SHong Zhang     /* Double the amount of total space in the list */
90fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
91fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
92fb908031SHong Zhang       ndouble++;
93fb908031SHong Zhang     }
94fb908031SHong Zhang 
95fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
96fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
97fb908031SHong Zhang     current_space->array           += cnzi;
98fb908031SHong Zhang     current_space->local_used      += cnzi;
99fb908031SHong Zhang     current_space->local_remaining -= cnzi;
100fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
101fb908031SHong Zhang   }
102fb908031SHong Zhang 
103fb908031SHong Zhang   /* Column indices are in the list of free space */
104fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
105fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
106fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
107fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
108fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
109b561aa9dSHong Zhang 
11026be0446SHong Zhang   /* put together the new symbolic matrix */
111aa1aec99SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr);
11258c24d83SHong Zhang 
11358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
11458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
11558c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
116aa1aec99SHong Zhang   c->free_a  = PETSC_FALSE;
117e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
11858c24d83SHong Zhang   c->nonew   = 0;
119002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1200b7e3e3dSHong Zhang 
1217212da7cSHong Zhang   /* set MatInfo */
1227212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
123f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1247212da7cSHong Zhang   c->maxnz                     = ci[am];
1257212da7cSHong Zhang   c->nz                        = ci[am];
126fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1277212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1287212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1297212da7cSHong Zhang 
1307212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1317212da7cSHong Zhang   if (ci[am]) {
132fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
133f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
134f2b054eeSHong Zhang   } else {
135f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
136be0fcf8dSHong Zhang   }
137f2b054eeSHong Zhang #endif
13858c24d83SHong Zhang   PetscFunctionReturn(0);
13958c24d83SHong Zhang }
140d50806bdSBarry Smith 
14126be0446SHong Zhang #undef __FUNCT__
14226be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
143dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
144d50806bdSBarry Smith {
145dfbe8321SBarry Smith   PetscErrorCode ierr;
146d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
147d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
148d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
149d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
15038baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
151c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
152519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
153aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
154aa1aec99SHong Zhang   PetscScalar    *ab_dense;
155d50806bdSBarry Smith 
156d50806bdSBarry Smith   PetscFunctionBegin;
157aa1aec99SHong Zhang   /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */
158aa1aec99SHong Zhang   if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
159aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
160aa1aec99SHong Zhang     c->a      = ca;
161aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
162aa1aec99SHong Zhang 
163aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
164aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
165aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
166aa1aec99SHong Zhang   } else {
167aa1aec99SHong Zhang     ca       = c->a;
168aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
169aa1aec99SHong Zhang   }
170aa1aec99SHong Zhang 
171c124e916SHong Zhang   /* clean old values in C */
172c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
173d50806bdSBarry Smith   /* Traverse A row-wise. */
174d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
175d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
176d50806bdSBarry Smith   for (i=0; i<am; i++) {
177d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
178d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
179519eb980SHong Zhang       brow = aj[j];
180d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
181d50806bdSBarry Smith       bjj  = bj + bi[brow];
182d50806bdSBarry Smith       baj  = ba + bi[brow];
183519eb980SHong Zhang       /* perform dense axpy */
18436ec6d2dSHong Zhang       valtmp = aa[j];
185519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
18636ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
187519eb980SHong Zhang       }
188519eb980SHong Zhang       flops += 2*bnzi;
189519eb980SHong Zhang     }
190c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
191c58ca2e3SHong Zhang 
192c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
193519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
194519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
195519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
196519eb980SHong Zhang     }
197519eb980SHong Zhang     flops += cnzi;
198519eb980SHong Zhang     cj += cnzi; ca += cnzi;
199519eb980SHong Zhang   }
200c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
201c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
202c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
203c58ca2e3SHong Zhang   PetscFunctionReturn(0);
204c58ca2e3SHong Zhang }
205c58ca2e3SHong Zhang 
206c58ca2e3SHong Zhang #undef __FUNCT__
20725023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
20825023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
209c58ca2e3SHong Zhang {
210c58ca2e3SHong Zhang   PetscErrorCode ierr;
211c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
212c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
213c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
214c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
215c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
216c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
217c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
21836ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
219c58ca2e3SHong Zhang   PetscInt       nextb;
220c58ca2e3SHong Zhang 
221c58ca2e3SHong Zhang   PetscFunctionBegin;
222e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
223971236abSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr);
224e9e4536cSHong Zhang #endif
225c58ca2e3SHong Zhang   /* clean old values in C */
226c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
227c58ca2e3SHong Zhang   /* Traverse A row-wise. */
228c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
229c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
230519eb980SHong Zhang   for (i=0;i<am;i++) {
231519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
232519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
233519eb980SHong Zhang     for (j=0;j<anzi;j++) {
234519eb980SHong Zhang       brow = aj[j];
235519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
236519eb980SHong Zhang       bjj  = bj + bi[brow];
237519eb980SHong Zhang       baj  = ba + bi[brow];
238519eb980SHong Zhang       /* perform sparse axpy */
23936ec6d2dSHong Zhang       valtmp = aa[j];
240c124e916SHong Zhang       nextb  = 0;
241c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
242c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
24336ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
244c124e916SHong Zhang         }
245d50806bdSBarry Smith       }
246d50806bdSBarry Smith       flops += 2*bnzi;
247d50806bdSBarry Smith     }
248519eb980SHong Zhang     aj += anzi; aa += anzi;
249519eb980SHong Zhang     cj += cnzi; ca += cnzi;
250d50806bdSBarry Smith   }
251c58ca2e3SHong Zhang 
252716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
255d50806bdSBarry Smith   PetscFunctionReturn(0);
256d50806bdSBarry Smith }
257bc011b1eSHong Zhang 
258e9e4536cSHong Zhang #undef __FUNCT__
2593c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2603c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
26125296bd5SBarry Smith {
26225296bd5SBarry Smith   PetscErrorCode     ierr;
26325296bd5SBarry Smith   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
26425296bd5SBarry Smith   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
2653c50cad2SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
26625296bd5SBarry Smith   MatScalar          *ca;
26725296bd5SBarry Smith   PetscReal          afill;
2683c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2693c50cad2SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
27025296bd5SBarry Smith 
27125296bd5SBarry Smith   PetscFunctionBegin;
2723c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2733c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2743c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2753c50cad2SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2763c50cad2SHong Zhang   ci[0] = 0;
2773c50cad2SHong Zhang 
2783c50cad2SHong Zhang   /* create and initialize a linked list */
2793c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2803c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2813c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
2823c50cad2SHong Zhang 
2833c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
2843c50cad2SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
2853c50cad2SHong Zhang   current_space = free_space;
2863c50cad2SHong Zhang 
2873c50cad2SHong Zhang   /* Determine ci and cj */
2883c50cad2SHong Zhang   for (i=0; i<am; i++) {
2893c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
2903c50cad2SHong Zhang     aj   = a->j + ai[i];
2913c50cad2SHong Zhang     for (j=0; j<anzi; j++){
2923c50cad2SHong Zhang       brow = aj[j];
2933c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
2943c50cad2SHong Zhang       bj   = b->j + bi[brow];
2953c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
2963c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
2973c50cad2SHong Zhang     }
2983c50cad2SHong Zhang     cnzi = lnk[1];
2993c50cad2SHong Zhang 
3003c50cad2SHong Zhang     /* If free space is not available, make more free space */
3013c50cad2SHong Zhang     /* Double the amount of total space in the list */
3023c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3033c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3043c50cad2SHong Zhang       ndouble++;
3053c50cad2SHong Zhang     }
3063c50cad2SHong Zhang 
3073c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3083c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3093c50cad2SHong Zhang     current_space->array           += cnzi;
3103c50cad2SHong Zhang     current_space->local_used      += cnzi;
3113c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3123c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3133c50cad2SHong Zhang   }
3143c50cad2SHong Zhang 
3153c50cad2SHong Zhang   /* Column indices are in the list of free space */
3163c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3173c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3183c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3193c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3203c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
32125296bd5SBarry Smith 
32225296bd5SBarry Smith   /* Allocate space for ca */
32325296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
32425296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
32525296bd5SBarry Smith 
32625296bd5SBarry Smith   /* put together the new symbolic matrix */
32725296bd5SBarry Smith   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
32825296bd5SBarry Smith 
32925296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
33025296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
33125296bd5SBarry Smith   c = (Mat_SeqAIJ *)((*C)->data);
33225296bd5SBarry Smith   c->free_a   = PETSC_TRUE;
33325296bd5SBarry Smith   c->free_ij  = PETSC_TRUE;
33425296bd5SBarry Smith   c->nonew    = 0;
33525296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
33625296bd5SBarry Smith 
33725296bd5SBarry Smith   /* set MatInfo */
33825296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
33925296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
34025296bd5SBarry Smith   c->maxnz                     = ci[am];
34125296bd5SBarry Smith   c->nz                        = ci[am];
3423c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
34325296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
34425296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
34525296bd5SBarry Smith 
34625296bd5SBarry Smith #if defined(PETSC_USE_INFO)
34725296bd5SBarry Smith   if (ci[am]) {
3483c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
34925296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
35025296bd5SBarry Smith   } else {
35125296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
35225296bd5SBarry Smith   }
35325296bd5SBarry Smith #endif
35425296bd5SBarry Smith   PetscFunctionReturn(0);
35525296bd5SBarry Smith }
35625296bd5SBarry Smith 
35725296bd5SBarry Smith 
35825296bd5SBarry Smith #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;
36525c41797SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
366e9e4536cSHong Zhang   MatScalar          *ca;
367e9e4536cSHong Zhang   PetscReal          afill;
368bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
369bd958071SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
370e9e4536cSHong Zhang 
371e9e4536cSHong Zhang   PetscFunctionBegin;
372bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
373bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
374bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
375bd958071SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
376bd958071SHong Zhang   ci[0] = 0;
377bd958071SHong Zhang 
378bd958071SHong Zhang   /* create and initialize a linked list */
379bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
380bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
381bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
382bd958071SHong Zhang 
383bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
384bd958071SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
385bd958071SHong Zhang   current_space = free_space;
386bd958071SHong Zhang 
387bd958071SHong Zhang   /* Determine ci and cj */
388bd958071SHong Zhang   for (i=0; i<am; i++) {
389bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
390bd958071SHong Zhang     aj   = a->j + ai[i];
391bd958071SHong Zhang     for (j=0; j<anzi; j++){
392bd958071SHong Zhang       brow = aj[j];
393bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
394bd958071SHong Zhang       bj   = b->j + bi[brow];
395bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
396bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
397bd958071SHong Zhang     }
398bd958071SHong Zhang     cnzi = lnk[0];
399bd958071SHong Zhang 
400bd958071SHong Zhang     /* If free space is not available, make more free space */
401bd958071SHong Zhang     /* Double the amount of total space in the list */
402bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
403bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
404bd958071SHong Zhang       ndouble++;
405bd958071SHong Zhang     }
406bd958071SHong Zhang 
407bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
408bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
409bd958071SHong Zhang     current_space->array           += cnzi;
410bd958071SHong Zhang     current_space->local_used      += cnzi;
411bd958071SHong Zhang     current_space->local_remaining -= cnzi;
412bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
413bd958071SHong Zhang   }
414bd958071SHong Zhang 
415bd958071SHong Zhang   /* Column indices are in the list of free space */
416bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
417bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
418bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
419bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
420bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
421e9e4536cSHong Zhang 
422e9e4536cSHong Zhang   /* Allocate space for ca */
423bd958071SHong Zhang   /*-----------------------*/
424e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
425e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
426e9e4536cSHong Zhang 
427e9e4536cSHong Zhang   /* put together the new symbolic matrix */
428e9e4536cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
429e9e4536cSHong Zhang 
430e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
431e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
432e9e4536cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
433e9e4536cSHong Zhang   c->free_a   = PETSC_TRUE;
434e9e4536cSHong Zhang   c->free_ij  = PETSC_TRUE;
435e9e4536cSHong Zhang   c->nonew    = 0;
43625023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
437e9e4536cSHong Zhang 
438e9e4536cSHong Zhang   /* set MatInfo */
439e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
440e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
441e9e4536cSHong Zhang   c->maxnz                     = ci[am];
442e9e4536cSHong Zhang   c->nz                        = ci[am];
443bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
444e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
445e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
446e9e4536cSHong Zhang 
447e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
448e9e4536cSHong Zhang   if (ci[am]) {
449bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
450e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
451e9e4536cSHong Zhang   } else {
452e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
453e9e4536cSHong Zhang   }
454e9e4536cSHong Zhang #endif
455e9e4536cSHong Zhang   PetscFunctionReturn(0);
456e9e4536cSHong Zhang }
457e9e4536cSHong Zhang 
458e9e4536cSHong Zhang 
459d2853540SHong Zhang /* This routine is not used. Should be removed! */
460bc011b1eSHong Zhang #undef __FUNCT__
4616fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
4626fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4635df89d91SHong Zhang {
464bc011b1eSHong Zhang   PetscErrorCode ierr;
465bc011b1eSHong Zhang 
466bc011b1eSHong Zhang   PetscFunctionBegin;
467bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
4686fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
469bc011b1eSHong Zhang   }
4706fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
471bc011b1eSHong Zhang   PetscFunctionReturn(0);
472bc011b1eSHong Zhang }
473bc011b1eSHong Zhang 
474bc011b1eSHong Zhang #undef __FUNCT__
475e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
476e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
4772128a86cSHong Zhang {
4782128a86cSHong Zhang   PetscErrorCode      ierr;
479e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
4802128a86cSHong Zhang 
4812128a86cSHong Zhang   PetscFunctionBegin;
482b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
4832128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
4842128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
4852128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
4862128a86cSHong Zhang   PetscFunctionReturn(0);
4872128a86cSHong Zhang }
4882128a86cSHong Zhang 
4892128a86cSHong Zhang #undef __FUNCT__
4902128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
4912128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
4922128a86cSHong Zhang {
4932128a86cSHong Zhang   PetscErrorCode      ierr;
4942128a86cSHong Zhang   PetscContainer      container;
495e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
4962128a86cSHong Zhang 
4972128a86cSHong Zhang   PetscFunctionBegin;
498e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
4992128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
5002128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
5012128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
5022128a86cSHong Zhang   if (A->ops->destroy) {
5032128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
5042128a86cSHong Zhang   }
505e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
5062128a86cSHong Zhang   PetscFunctionReturn(0);
5072128a86cSHong Zhang }
5082128a86cSHong Zhang 
5092128a86cSHong Zhang #undef __FUNCT__
5106fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
5116fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
512bc011b1eSHong Zhang {
5135df89d91SHong Zhang   PetscErrorCode      ierr;
51481d82fe4SHong Zhang   Mat                 Bt;
51581d82fe4SHong Zhang   PetscInt            *bti,*btj;
516e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
5172128a86cSHong Zhang   PetscContainer      container;
518d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
519d2853540SHong Zhang 
52081d82fe4SHong Zhang   PetscFunctionBegin;
521d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
52281d82fe4SHong Zhang    /* create symbolic Bt */
52381d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
52481d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
52581d82fe4SHong Zhang 
52681d82fe4SHong Zhang   /* get symbolic C=A*Bt */
52781d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
52881d82fe4SHong Zhang 
5292128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
530e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
5312128a86cSHong Zhang 
5322128a86cSHong Zhang   /* attach the supporting struct to C */
5332128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
5342128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
535e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
536e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
5372128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5382128a86cSHong Zhang 
5392128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
5402128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
5412128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
5422128a86cSHong Zhang 
543d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
544d2853540SHong Zhang   etime2 += tf - t0;
545d2853540SHong Zhang 
546f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
5472128a86cSHong Zhang   if (multtrans->usecoloring){
548b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
549b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
5502128a86cSHong Zhang     ISColoring           iscoloring;
5512128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
552d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
553e8354b3eSHong Zhang 
554e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
555502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
556d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
557d2853540SHong Zhang     etime0 += tf - t0;
558d2853540SHong Zhang 
559d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
560b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
5612128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
5622128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
563e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
564d2853540SHong Zhang     etime01 += tf - t0;
5652128a86cSHong Zhang 
566e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
5672128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
5682128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
5692128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
5702128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
5712128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
5722128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
5732128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
5742128a86cSHong Zhang 
5752128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
5762128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
5772128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
5782128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
5792128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
5802128a86cSHong Zhang     multtrans->ABt_den = C_dense;
581e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
582e8354b3eSHong Zhang     etime1 += tf - t0;
583f94ccd6cSHong Zhang 
584f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
5851ce71dffSSatish Balay     {
586f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
587f94ccd6cSHong 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));
588f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
5891ce71dffSSatish Balay     }
590f94ccd6cSHong Zhang #endif
5912128a86cSHong Zhang   }
592e99be685SHong Zhang   /* clean up */
593e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
594e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
5952128a86cSHong Zhang 
596f94ccd6cSHong Zhang 
597f94ccd6cSHong Zhang 
59881d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
59981d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
6001fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
6011fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6021fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
6031fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
6041fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
6051fa1209cSHong Zhang   MatScalar          *ca;
6061fa1209cSHong Zhang   PetscBT            lnkbt;
6071fa1209cSHong Zhang   PetscReal          afill;
6085df89d91SHong Zhang 
6091fa1209cSHong Zhang   /* Allocate row pointer array ci  */
6101fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6111fa1209cSHong Zhang   ci[0] = 0;
6121fa1209cSHong Zhang 
6131fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
6141fa1209cSHong Zhang   nlnk = bm+1;
6151fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
6161fa1209cSHong Zhang 
6171fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
6181fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6191fa1209cSHong Zhang   current_space = free_space;
6201fa1209cSHong Zhang 
6211fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
6221fa1209cSHong Zhang   for (i=0; i<am; i++) {
6231fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
6241fa1209cSHong Zhang     cnzi = 0;
6251fa1209cSHong Zhang     acol = aj + ai[i];
6261fa1209cSHong Zhang     for (j=0; j<bm; j++){
6271fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
6281fa1209cSHong Zhang       bcol= bj + bi[j];
62981d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
6301fa1209cSHong Zhang       ka = 0; kb = 0;
6311fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
63281d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
63381d82fe4SHong Zhang         if (ka == anzi) break;
63481d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
63581d82fe4SHong Zhang         if (kb == bnzj) break;
6361fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
6371fa1209cSHong Zhang           index[0] = j;
6381fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
6391fa1209cSHong Zhang           cnzi++;
6401fa1209cSHong Zhang           break;
6411fa1209cSHong Zhang         }
6421fa1209cSHong Zhang       }
6431fa1209cSHong Zhang     }
6441fa1209cSHong Zhang 
6451fa1209cSHong Zhang     /* If free space is not available, make more free space */
6461fa1209cSHong Zhang     /* Double the amount of total space in the list */
6471fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
6481fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
6491fa1209cSHong Zhang       nspacedouble++;
6501fa1209cSHong Zhang     }
6511fa1209cSHong Zhang 
6521fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
6531fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
6541fa1209cSHong Zhang     current_space->array           += cnzi;
6551fa1209cSHong Zhang     current_space->local_used      += cnzi;
6561fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
6571fa1209cSHong Zhang 
6581fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
6591fa1209cSHong Zhang   }
6601fa1209cSHong Zhang 
6611fa1209cSHong Zhang 
6621fa1209cSHong Zhang   /* Column indices are in the list of free space.
6631fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
6641fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6651fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6661fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
6671fa1209cSHong Zhang 
6681fa1209cSHong Zhang   /* Allocate space for ca */
6691fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
6701fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
6711fa1209cSHong Zhang 
6721fa1209cSHong Zhang   /* put together the new symbolic matrix */
6731fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
6741fa1209cSHong Zhang 
6751fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6761fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6771fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
6781fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
6791fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
6801fa1209cSHong Zhang   c->nonew    = 0;
6811fa1209cSHong Zhang 
6821fa1209cSHong Zhang   /* set MatInfo */
6831fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6841fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
6851fa1209cSHong Zhang   c->maxnz                     = ci[am];
6861fa1209cSHong Zhang   c->nz                        = ci[am];
6871fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
6881fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
6891fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
6901fa1209cSHong Zhang 
6911fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
6921fa1209cSHong Zhang   if (ci[am]) {
6931fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
6946fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
6951fa1209cSHong Zhang   } else {
6961fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6971fa1209cSHong Zhang   }
6981fa1209cSHong Zhang #endif
69981d82fe4SHong Zhang #endif
7005df89d91SHong Zhang   PetscFunctionReturn(0);
7015df89d91SHong Zhang }
7025df89d91SHong Zhang 
7036973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
7045df89d91SHong Zhang #undef __FUNCT__
7056fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
7066fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
7075df89d91SHong Zhang {
7085df89d91SHong Zhang   PetscErrorCode ierr;
7095df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
710e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
71181d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
7125df89d91SHong Zhang   PetscLogDouble flops=0.0;
713aa1aec99SHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval;
714e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
7152128a86cSHong Zhang   PetscContainer      container;
7166973769bSHong Zhang #if defined(USE_ARRAY)
7176973769bSHong Zhang   MatScalar      *spdot;
7186973769bSHong Zhang #endif
7195df89d91SHong Zhang 
7205df89d91SHong Zhang   PetscFunctionBegin;
721e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
7222128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
7232128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
7242128a86cSHong Zhang   if (multtrans->usecoloring){
725b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
726c8db22f6SHong Zhang     Mat                   Bt_dense;
727c8db22f6SHong Zhang     PetscInt              m,n;
728b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
729a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
730a0b95ffbSSatish Balay 
7312128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
732c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
733c8db22f6SHong Zhang 
734b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
7352128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
736b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
7372128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
738b2d2b04fSHong Zhang     etime0 += tf - t0;
739c8db22f6SHong Zhang 
740c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
7412128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
7422128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
7432128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7442128a86cSHong Zhang     etime2 += tf - t0;
745c8db22f6SHong Zhang 
7462128a86cSHong Zhang     /* Recover C from C_dense */
7472128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
748b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
7492128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7502128a86cSHong Zhang     etime1 += tf - t0;
751f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
752f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
753f94ccd6cSHong Zhang #endif
754c8db22f6SHong Zhang     PetscFunctionReturn(0);
755c8db22f6SHong Zhang   }
756c8db22f6SHong Zhang 
7576973769bSHong Zhang #if defined(USE_ARRAY)
7586973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
759e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
760e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
7616973769bSHong Zhang #endif
7626973769bSHong Zhang 
76381d82fe4SHong Zhang   /* clear old values in C */
764aa1aec99SHong Zhang   if (!c->a){
765aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
766aa1aec99SHong Zhang     c->a      = ca;
767aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
768aa1aec99SHong Zhang   } else {
769aa1aec99SHong Zhang     ca =  c->a;
770aa1aec99SHong Zhang   }
77181d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
7725df89d91SHong Zhang 
7731fa1209cSHong Zhang   for (i=0; i<cm; i++) {
77481d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
77581d82fe4SHong Zhang     acol = aj + ai[i];
7766973769bSHong Zhang     aval = aa + ai[i];
7771fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
7781fa1209cSHong Zhang     ccol = cj + ci[i];
7796973769bSHong Zhang     cval = ca + ci[i];
7801fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
78181d82fe4SHong Zhang       brow = ccol[j];
78281d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
78381d82fe4SHong Zhang       bcol = bj + bi[brow];
7846973769bSHong Zhang       bval = ba + bi[brow];
7856973769bSHong Zhang 
78681d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
7876973769bSHong Zhang #if defined(USE_ARRAY)
7886973769bSHong Zhang       /* put ba to spdot array */
7896973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
7906973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
7916973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
7926973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
7936973769bSHong Zhang       }
7946973769bSHong Zhang       /* zero spdot array */
7956973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
7966973769bSHong Zhang #else
79781d82fe4SHong Zhang       nexta = 0; nextb = 0;
79881d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
79981d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
80081d82fe4SHong Zhang         if (nexta == anzi) break;
80181d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
80281d82fe4SHong Zhang         if (nextb == bnzj) break;
80381d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
8046973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
80581d82fe4SHong Zhang           nexta++; nextb++;
80681d82fe4SHong Zhang           flops += 2;
8071fa1209cSHong Zhang         }
8081fa1209cSHong Zhang       }
8096973769bSHong Zhang #endif
81081d82fe4SHong Zhang     }
81181d82fe4SHong Zhang   }
81281d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81381d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81481d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
8156973769bSHong Zhang #if defined(USE_ARRAY)
8166973769bSHong Zhang   ierr = PetscFree(spdot);
8176973769bSHong Zhang #endif
8185df89d91SHong Zhang   PetscFunctionReturn(0);
8195df89d91SHong Zhang }
8205df89d91SHong Zhang 
8215df89d91SHong Zhang #undef __FUNCT__
82275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
82375648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
8245df89d91SHong Zhang   PetscErrorCode ierr;
8255df89d91SHong Zhang 
8265df89d91SHong Zhang   PetscFunctionBegin;
8275df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
82875648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8295df89d91SHong Zhang   }
83075648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8315df89d91SHong Zhang   PetscFunctionReturn(0);
8325df89d91SHong Zhang }
8335df89d91SHong Zhang 
8345df89d91SHong Zhang #undef __FUNCT__
83575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
83675648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
8375df89d91SHong Zhang {
838bc011b1eSHong Zhang   PetscErrorCode ierr;
839bc011b1eSHong Zhang   Mat            At;
84038baddfdSBarry Smith   PetscInt       *ati,*atj;
841bc011b1eSHong Zhang 
842bc011b1eSHong Zhang   PetscFunctionBegin;
843bc011b1eSHong Zhang   /* create symbolic At */
844bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
845d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
846bc011b1eSHong Zhang 
847bc011b1eSHong Zhang   /* get symbolic C=At*B */
848bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
849bc011b1eSHong Zhang 
850bc011b1eSHong Zhang   /* clean up */
8516bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
852bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
853bc011b1eSHong Zhang   PetscFunctionReturn(0);
854bc011b1eSHong Zhang }
855bc011b1eSHong Zhang 
856bc011b1eSHong Zhang #undef __FUNCT__
85775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
85875648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
859bc011b1eSHong Zhang {
860bc011b1eSHong Zhang   PetscErrorCode ierr;
8610fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
862d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
863d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
864d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
865aa1aec99SHong Zhang   MatScalar      *aa=a->a,*ba,*ca,*caj;
866bc011b1eSHong Zhang 
867bc011b1eSHong Zhang   PetscFunctionBegin;
868aa1aec99SHong Zhang   if (!c->a){
869aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
870aa1aec99SHong Zhang     c->a      = ca;
871aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
872aa1aec99SHong Zhang   } else {
873aa1aec99SHong Zhang     ca = c->a;
874aa1aec99SHong Zhang   }
875bc011b1eSHong Zhang   /* clear old values in C */
876bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
877bc011b1eSHong Zhang 
878bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
879bc011b1eSHong Zhang   for (i=0;i<am;i++) {
880bc011b1eSHong Zhang     bj   = b->j + bi[i];
881bc011b1eSHong Zhang     ba   = b->a + bi[i];
882bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
883bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
884bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
885bc011b1eSHong Zhang       nextb = 0;
8860fbc74f4SHong Zhang       crow  = *aj++;
887bc011b1eSHong Zhang       cjj   = cj + ci[crow];
888bc011b1eSHong Zhang       caj   = ca + ci[crow];
889bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
890bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
8910fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
8920fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
893bc011b1eSHong Zhang           nextb++;
894bc011b1eSHong Zhang         }
895bc011b1eSHong Zhang       }
896bc011b1eSHong Zhang       flops += 2*bnzi;
8970fbc74f4SHong Zhang       aa++;
898bc011b1eSHong Zhang     }
899bc011b1eSHong Zhang   }
900bc011b1eSHong Zhang 
901bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
902bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
903bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
904bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
905bc011b1eSHong Zhang   PetscFunctionReturn(0);
906bc011b1eSHong Zhang }
9079123193aSHong Zhang 
908ec01deb9SMatthew Knepley EXTERN_C_BEGIN
9099123193aSHong Zhang #undef __FUNCT__
9109123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
9119123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9129123193aSHong Zhang {
9139123193aSHong Zhang   PetscErrorCode ierr;
9149123193aSHong Zhang 
9159123193aSHong Zhang   PetscFunctionBegin;
9169123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
9179123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
9189123193aSHong Zhang   }
9199123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
9209123193aSHong Zhang   PetscFunctionReturn(0);
9219123193aSHong Zhang }
922ec01deb9SMatthew Knepley EXTERN_C_END
923edd81eecSMatthew Knepley 
9249123193aSHong Zhang #undef __FUNCT__
9259123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
9269123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
9279123193aSHong Zhang {
9289123193aSHong Zhang   PetscErrorCode ierr;
9299123193aSHong Zhang 
9309123193aSHong Zhang   PetscFunctionBegin;
9315a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
9328cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
9339123193aSHong Zhang   PetscFunctionReturn(0);
9349123193aSHong Zhang }
9359123193aSHong Zhang 
9369123193aSHong Zhang #undef __FUNCT__
9379123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
9389123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
9399123193aSHong Zhang {
940f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
9419123193aSHong Zhang   PetscErrorCode ierr;
942dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
943dd6ea824SBarry Smith   MatScalar      *aa;
944d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
945f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
9469123193aSHong Zhang 
9479123193aSHong Zhang   PetscFunctionBegin;
948f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
949e32f2f54SBarry 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);
950e32f2f54SBarry 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);
951e32f2f54SBarry 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);
952f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
953f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
954f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
955f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
956f32d5d43SBarry Smith     colam = col*am;
957f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
958f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
959f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
960f32d5d43SBarry Smith       aj  = a->j + a->i[i];
961f32d5d43SBarry Smith       aa  = a->a + a->i[i];
962f32d5d43SBarry Smith       for (j=0; j<n; j++) {
963f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
964f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
965f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
966f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
9679123193aSHong Zhang       }
968f32d5d43SBarry Smith       c[colam + i]       = r1;
969f32d5d43SBarry Smith       c[colam + am + i]  = r2;
970f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
971f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
972f32d5d43SBarry Smith     }
973f32d5d43SBarry Smith     b1 += bm4;
974f32d5d43SBarry Smith     b2 += bm4;
975f32d5d43SBarry Smith     b3 += bm4;
976f32d5d43SBarry Smith     b4 += bm4;
977f32d5d43SBarry Smith   }
978f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
979f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
980f32d5d43SBarry Smith       r1 = 0.0;
981f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
982f32d5d43SBarry Smith       aj  = a->j + a->i[i];
983f32d5d43SBarry Smith       aa  = a->a + a->i[i];
984f32d5d43SBarry Smith 
985f32d5d43SBarry Smith       for (j=0; j<n; j++) {
986f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
987f32d5d43SBarry Smith       }
988f32d5d43SBarry Smith       c[col*am + i]     = r1;
989f32d5d43SBarry Smith     }
990f32d5d43SBarry Smith     b1 += bm;
991f32d5d43SBarry Smith   }
992dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
993f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
994f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
995f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
996f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
997f32d5d43SBarry Smith   PetscFunctionReturn(0);
998f32d5d43SBarry Smith }
999f32d5d43SBarry Smith 
1000f32d5d43SBarry Smith /*
10014324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1002f32d5d43SBarry Smith */
1003f32d5d43SBarry Smith #undef __FUNCT__
1004f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1005f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1006f32d5d43SBarry Smith {
1007f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1008f32d5d43SBarry Smith   PetscErrorCode ierr;
1009dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1010dd6ea824SBarry Smith   MatScalar      *aa;
1011d0f46423SBarry 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;
10124324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1013f32d5d43SBarry Smith 
1014f32d5d43SBarry Smith   PetscFunctionBegin;
1015f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1016f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
1017f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
1018f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
10194324174eSBarry Smith 
10204324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
10214324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
10224324174eSBarry Smith       colam = col*am;
10234324174eSBarry Smith       arm   = a->compressedrow.nrows;
10244324174eSBarry Smith       ii    = a->compressedrow.i;
10254324174eSBarry Smith       ridx  = a->compressedrow.rindex;
10264324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
10274324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
10284324174eSBarry Smith 	n   = ii[i+1] - ii[i];
10294324174eSBarry Smith 	aj  = a->j + ii[i];
10304324174eSBarry Smith 	aa  = a->a + ii[i];
10314324174eSBarry Smith 	for (j=0; j<n; j++) {
10324324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
10334324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
10344324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
10354324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
10364324174eSBarry Smith 	}
10374324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
10384324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
10394324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
10404324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
10414324174eSBarry Smith       }
10424324174eSBarry Smith       b1 += bm4;
10434324174eSBarry Smith       b2 += bm4;
10444324174eSBarry Smith       b3 += bm4;
10454324174eSBarry Smith       b4 += bm4;
10464324174eSBarry Smith     }
10474324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
10484324174eSBarry Smith       colam = col*am;
10494324174eSBarry Smith       arm   = a->compressedrow.nrows;
10504324174eSBarry Smith       ii    = a->compressedrow.i;
10514324174eSBarry Smith       ridx  = a->compressedrow.rindex;
10524324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
10534324174eSBarry Smith 	r1 = 0.0;
10544324174eSBarry Smith 	n   = ii[i+1] - ii[i];
10554324174eSBarry Smith 	aj  = a->j + ii[i];
10564324174eSBarry Smith 	aa  = a->a + ii[i];
10574324174eSBarry Smith 
10584324174eSBarry Smith 	for (j=0; j<n; j++) {
10594324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
10604324174eSBarry Smith 	}
10614324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
10624324174eSBarry Smith       }
10634324174eSBarry Smith       b1 += bm;
10644324174eSBarry Smith     }
10654324174eSBarry Smith   } else {
1066f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
1067f32d5d43SBarry Smith       colam = col*am;
1068f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1069f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
1070f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1071f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1072f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1073f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1074f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
1075f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
1076f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
1077f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
1078f32d5d43SBarry Smith 	}
1079f32d5d43SBarry Smith 	c[colam + i]       += r1;
1080f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
1081f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
1082f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
1083f32d5d43SBarry Smith       }
1084f32d5d43SBarry Smith       b1 += bm4;
1085f32d5d43SBarry Smith       b2 += bm4;
1086f32d5d43SBarry Smith       b3 += bm4;
1087f32d5d43SBarry Smith       b4 += bm4;
1088f32d5d43SBarry Smith     }
1089f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
1090f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1091f32d5d43SBarry Smith 	r1 = 0.0;
1092f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1093f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1094f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1095f32d5d43SBarry Smith 
1096f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1097f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
1098f32d5d43SBarry Smith 	}
1099f32d5d43SBarry Smith 	c[col*am + i]     += r1;
1100f32d5d43SBarry Smith       }
1101f32d5d43SBarry Smith       b1 += bm;
1102f32d5d43SBarry Smith     }
11034324174eSBarry Smith   }
1104dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
1105f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1106f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
11079123193aSHong Zhang   PetscFunctionReturn(0);
11089123193aSHong Zhang }
1109b1683b59SHong Zhang 
1110b1683b59SHong Zhang #undef __FUNCT__
1111b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1112b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1113c8db22f6SHong Zhang {
1114c8db22f6SHong Zhang   PetscErrorCode ierr;
11152f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
11162f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
11172f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
11182f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
11192f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
11202f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1121c8db22f6SHong Zhang 
1122c8db22f6SHong Zhang   PetscFunctionBegin;
11232f5f1f90SHong Zhang   btval_den=btdense->v;
11242f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
11252f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
11262f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
11272f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1128d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
11292f5f1f90SHong Zhang       btcol = bj + bi[col];
11302f5f1f90SHong Zhang       btval = ba + bi[col];
11312f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1132d2853540SHong Zhang       for (j=0; j<anz; j++){
11332f5f1f90SHong Zhang         brow            = btcol[j];
11342f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1135c8db22f6SHong Zhang       }
1136c8db22f6SHong Zhang     }
11372f5f1f90SHong Zhang     btval_den += m;
1138c8db22f6SHong Zhang   }
1139c8db22f6SHong Zhang   PetscFunctionReturn(0);
1140c8db22f6SHong Zhang }
1141c8db22f6SHong Zhang 
1142c8db22f6SHong Zhang #undef __FUNCT__
1143b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1144b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
11458972f759SHong Zhang {
11468972f759SHong Zhang   PetscErrorCode ierr;
1147b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
11482f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1149b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
11502f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
11518972f759SHong Zhang 
11528972f759SHong Zhang   PetscFunctionBegin;
11538972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
1154b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
1155b2d2b04fSHong Zhang   cp_den = ca_den;
11562f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
11572f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
11582f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
11592f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
11602f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
11612f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1162b2d2b04fSHong Zhang     }
1163b2d2b04fSHong Zhang     cp_den += m;
1164b2d2b04fSHong Zhang   }
1165b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
11668972f759SHong Zhang   PetscFunctionReturn(0);
11678972f759SHong Zhang }
11688972f759SHong Zhang 
1169e99be685SHong Zhang /*
1170e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1171e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1172e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1173e99be685SHong Zhang  */
1174e99be685SHong Zhang #undef __FUNCT__
1175e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
1176e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1177e99be685SHong Zhang {
1178e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1179e99be685SHong Zhang   PetscErrorCode ierr;
1180e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1181e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1182e99be685SHong Zhang   PetscInt       *cspidx;
1183e99be685SHong Zhang 
1184e99be685SHong Zhang   PetscFunctionBegin;
1185e99be685SHong Zhang   *nn = n;
1186e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1187e99be685SHong Zhang   if (symmetric) {
1188e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
1189e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
1190e99be685SHong Zhang   } else {
1191e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1192e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1193e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1194e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1195e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1196e99be685SHong Zhang     jj = a->j;
1197e99be685SHong Zhang     for (i=0; i<nz; i++) {
1198e99be685SHong Zhang       collengths[jj[i]]++;
1199e99be685SHong Zhang     }
1200e99be685SHong Zhang     cia[0] = oshift;
1201e99be685SHong Zhang     for (i=0; i<n; i++) {
1202e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1203e99be685SHong Zhang     }
1204e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1205e99be685SHong Zhang     jj   = a->j;
1206e99be685SHong Zhang     for (row=0; row<m; row++) {
1207e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1208e99be685SHong Zhang       for (i=0; i<mr; i++) {
1209e99be685SHong Zhang         col = *jj++;
1210e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1211e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1212e99be685SHong Zhang       }
1213e99be685SHong Zhang     }
1214e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1215e99be685SHong Zhang     *ia = cia; *ja = cja;
1216e99be685SHong Zhang     *spidx = cspidx;
1217e99be685SHong Zhang   }
1218e99be685SHong Zhang   PetscFunctionReturn(0);
1219e99be685SHong Zhang }
1220e99be685SHong Zhang 
1221e99be685SHong Zhang #undef __FUNCT__
1222e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1223e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1224e99be685SHong Zhang {
1225e99be685SHong Zhang   PetscErrorCode ierr;
1226e99be685SHong Zhang 
1227e99be685SHong Zhang   PetscFunctionBegin;
1228e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1229e99be685SHong Zhang 
1230e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1231e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1232e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1233e99be685SHong Zhang   PetscFunctionReturn(0);
1234e99be685SHong Zhang }
1235e99be685SHong Zhang 
12368972f759SHong Zhang #undef __FUNCT__
1237b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1238b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1239b1683b59SHong Zhang {
1240b1683b59SHong Zhang   PetscErrorCode ierr;
1241d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1242b1683b59SHong Zhang   const PetscInt *is;
1243b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1244b1683b59SHong Zhang   IS             *isa;
1245b1683b59SHong Zhang   PetscBool      done;
1246b1683b59SHong Zhang   PetscBool      flg1,flg2;
1247b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1248e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1249d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1250b1683b59SHong Zhang 
1251b1683b59SHong Zhang   PetscFunctionBegin;
1252b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1253e99be685SHong Zhang 
1254b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1255*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1256*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1257b1683b59SHong Zhang   if (flg1 || flg2) {
1258b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1259b1683b59SHong Zhang   }
1260b1683b59SHong Zhang 
1261b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1262b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1263b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1264b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1265b1683b59SHong Zhang   c->rstart = 0;
1266b1683b59SHong Zhang 
1267b1683b59SHong Zhang   c->ncolors = nis;
1268b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1269b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1270d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1271d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1272d2853540SHong Zhang   colorforrow[0]    = 0;
1273d2853540SHong Zhang   rows_i            = rows;
1274d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1275b1683b59SHong Zhang 
1276d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1277d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1278d2853540SHong Zhang   colorforcol[0] = 0;
1279d2853540SHong Zhang   columns_i      = columns;
1280d2853540SHong Zhang 
1281e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1282b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1283b1683b59SHong Zhang 
1284b28d80bdSHong Zhang   cm = c->m;
1285b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1286e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1287b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1288b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1289b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1290b1683b59SHong Zhang     c->ncolumns[i] = n;
1291b1683b59SHong Zhang     if (n) {
1292d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1293b1683b59SHong Zhang     }
1294d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1295d2853540SHong Zhang     columns_i       += n;
1296b1683b59SHong Zhang 
1297b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1298e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1299e99be685SHong Zhang 
1300b1683b59SHong Zhang     /* loop over columns*/
1301b1683b59SHong Zhang     for (j=0; j<n; j++) {
1302b1683b59SHong Zhang       col     = is[j];
1303d2853540SHong Zhang       row_idx = cj + ci[col];
1304b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1305b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1306b1683b59SHong Zhang       for (k=0; k<m; k++) {
1307e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1308d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1309b1683b59SHong Zhang       }
1310b1683b59SHong Zhang     }
1311b1683b59SHong Zhang     /* count the number of hits */
1312b1683b59SHong Zhang     nrows = 0;
1313e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1314b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1315b1683b59SHong Zhang     }
1316b1683b59SHong Zhang     c->nrows[i]      = nrows;
1317d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1318d2853540SHong Zhang 
1319b1683b59SHong Zhang     nrows       = 0;
1320e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1321b1683b59SHong Zhang       if (rowhit[j]) {
1322d2853540SHong Zhang         rows_i[nrows]            = j;
1323e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1324b1683b59SHong Zhang         nrows++;
1325b1683b59SHong Zhang       }
1326b1683b59SHong Zhang     }
1327b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1328d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1329b1683b59SHong Zhang   }
1330e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1331b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1332b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1333d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1334d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1335d2853540SHong Zhang #endif
1336b28d80bdSHong Zhang 
1337d2853540SHong Zhang   c->colorforrow     = colorforrow;
1338d2853540SHong Zhang   c->rows            = rows;
1339d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1340d2853540SHong Zhang   c->colorforcol     = colorforcol;
1341d2853540SHong Zhang   c->columns         = columns;
1342f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1343b1683b59SHong Zhang   PetscFunctionReturn(0);
1344b1683b59SHong Zhang }
1345