xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision c373ccc68bfa07b9af08e54bc02c781257ff10f0)
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>
90ced3a2bSJed Brown #include <../src/mat/utils/petscheap.h>
10c6db04a5SJed Brown #include <petscbt.h>
11af0996ceSBarry Smith #include <petsc/private/isimpl.h>
1207475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h>
137bab7c10SHong Zhang 
1458cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
15cd093f1dSJed Brown 
166284ec50SHong Zhang #undef __FUNCT__
175c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
18150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1938baddfdSBarry Smith {
20dfbe8321SBarry Smith   PetscErrorCode ierr;
216540a9cdSHong Zhang   const char     *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
226540a9cdSHong Zhang   PetscInt       alg=0; /* set default algorithm */
235c66b693SKris Buschelman 
245c66b693SKris Buschelman   PetscFunctionBegin;
2526be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
26152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
276540a9cdSHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,6,algTypes[0],&alg,NULL);CHKERRQ(ierr);
28d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
293ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
306540a9cdSHong Zhang     switch (alg) {
316540a9cdSHong Zhang     case 1:
32aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
336540a9cdSHong Zhang       break;
346540a9cdSHong Zhang     case 2:
356540a9cdSHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
366540a9cdSHong Zhang       break;
376540a9cdSHong Zhang     case 3:
380ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
396540a9cdSHong Zhang       break;
406540a9cdSHong Zhang     case 4:
418a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
426540a9cdSHong Zhang       break;
436540a9cdSHong Zhang     case 5:
4458cf0668SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
456540a9cdSHong Zhang       break;
466540a9cdSHong Zhang     default:
4726be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
486540a9cdSHong Zhang      break;
4925023028SHong Zhang     }
503ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
5126be0446SHong Zhang   }
525c913ed7SHong Zhang 
533ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
5401e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
553ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
565c66b693SKris Buschelman   PetscFunctionReturn(0);
575c66b693SKris Buschelman }
581c24bd37SHong Zhang 
59e9e4536cSHong Zhang #undef __FUNCT__
6058cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed"
6158cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
62b561aa9dSHong Zhang {
63b561aa9dSHong Zhang   PetscErrorCode     ierr;
64b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
65971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
66c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
67b561aa9dSHong Zhang   PetscReal          afill;
68eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
69eca6b86aSHong Zhang   PetscTable         ta;
70fb908031SHong Zhang   PetscBT            lnkbt;
710298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
72b561aa9dSHong Zhang 
73b561aa9dSHong Zhang   PetscFunctionBegin;
74bd958071SHong Zhang   /* Get ci and cj */
75bd958071SHong Zhang   /*---------------*/
76fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
77fb908031SHong Zhang   /* free space for accumulating nonzero column info */
78854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
79fb908031SHong Zhang   ci[0] = 0;
80fb908031SHong Zhang 
81fb908031SHong Zhang   /* create and initialize a linked list */
82*c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
83*c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
84eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
85eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
86eca6b86aSHong Zhang 
87eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr);
88fb908031SHong Zhang 
89fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
90f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
912205254eSKarl Rupp 
92fb908031SHong Zhang   current_space = free_space;
93fb908031SHong Zhang 
94fb908031SHong Zhang   /* Determine ci and cj */
95fb908031SHong Zhang   for (i=0; i<am; i++) {
96fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
97fb908031SHong Zhang     aj   = a->j + ai[i];
98fb908031SHong Zhang     for (j=0; j<anzi; j++) {
99fb908031SHong Zhang       brow = aj[j];
100fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
101fb908031SHong Zhang       bj   = b->j + bi[brow];
102fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
103fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
104fb908031SHong Zhang     }
105fb908031SHong Zhang     cnzi = lnk[0];
106fb908031SHong Zhang 
107fb908031SHong Zhang     /* If free space is not available, make more free space */
108fb908031SHong Zhang     /* Double the amount of total space in the list */
109fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
110f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
111fb908031SHong Zhang       ndouble++;
112fb908031SHong Zhang     }
113fb908031SHong Zhang 
114fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
115fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1162205254eSKarl Rupp 
117fb908031SHong Zhang     current_space->array           += cnzi;
118fb908031SHong Zhang     current_space->local_used      += cnzi;
119fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1202205254eSKarl Rupp 
121fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
122fb908031SHong Zhang   }
123fb908031SHong Zhang 
124fb908031SHong Zhang   /* Column indices are in the list of free space */
125fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
126fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
127854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
128fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
129fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
130b561aa9dSHong Zhang 
13126be0446SHong Zhang   /* put together the new symbolic matrix */
132ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
13333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
13458c24d83SHong Zhang 
13558c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
13658c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
13758c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
138aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
139e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
14058c24d83SHong Zhang   c->nonew                  = 0;
141002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1420b7e3e3dSHong Zhang 
1437212da7cSHong Zhang   /* set MatInfo */
1447212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
145f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1467212da7cSHong Zhang   c->maxnz                     = ci[am];
1477212da7cSHong Zhang   c->nz                        = ci[am];
148fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1497212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1507212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1517212da7cSHong Zhang 
1527212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1537212da7cSHong Zhang   if (ci[am]) {
15457622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
15557622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
156f2b054eeSHong Zhang   } else {
157f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
158be0fcf8dSHong Zhang   }
159f2b054eeSHong Zhang #endif
16058c24d83SHong Zhang   PetscFunctionReturn(0);
16158c24d83SHong Zhang }
162d50806bdSBarry Smith 
16326be0446SHong Zhang #undef __FUNCT__
16426be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
165dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
166d50806bdSBarry Smith {
167dfbe8321SBarry Smith   PetscErrorCode ierr;
168d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
169d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
170d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
171d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
17238baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
173c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
174519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
175aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
176aa1aec99SHong Zhang   PetscScalar    *ab_dense;
177d50806bdSBarry Smith 
178d50806bdSBarry Smith   PetscFunctionBegin;
179aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
180854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
181aa1aec99SHong Zhang     c->a      = ca;
182aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
183aa1aec99SHong Zhang   } else {
184aa1aec99SHong Zhang     ca        = c->a;
185d90d86d0SMatthew G. Knepley   }
186d90d86d0SMatthew G. Knepley   if (!c->matmult_abdense) {
1871795a4d1SJed Brown     ierr = PetscCalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr);
188d90d86d0SMatthew G. Knepley     c->matmult_abdense = ab_dense;
189d90d86d0SMatthew G. Knepley   } else {
190aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
191aa1aec99SHong Zhang   }
192aa1aec99SHong Zhang 
193c124e916SHong Zhang   /* clean old values in C */
194c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
195d50806bdSBarry Smith   /* Traverse A row-wise. */
196d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
197d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
198d50806bdSBarry Smith   for (i=0; i<am; i++) {
199d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
200d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
201519eb980SHong Zhang       brow = aj[j];
202d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
203d50806bdSBarry Smith       bjj  = bj + bi[brow];
204d50806bdSBarry Smith       baj  = ba + bi[brow];
205519eb980SHong Zhang       /* perform dense axpy */
20636ec6d2dSHong Zhang       valtmp = aa[j];
207519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
20836ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
209519eb980SHong Zhang       }
210519eb980SHong Zhang       flops += 2*bnzi;
211519eb980SHong Zhang     }
212c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
213c58ca2e3SHong Zhang 
214c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
215519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
216519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
217519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
218519eb980SHong Zhang     }
219519eb980SHong Zhang     flops += cnzi;
220519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
221519eb980SHong Zhang   }
222c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
224c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
225c58ca2e3SHong Zhang   PetscFunctionReturn(0);
226c58ca2e3SHong Zhang }
227c58ca2e3SHong Zhang 
228c58ca2e3SHong Zhang #undef __FUNCT__
22925023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
23025023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
231c58ca2e3SHong Zhang {
232c58ca2e3SHong Zhang   PetscErrorCode ierr;
233c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
234c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
235c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
236c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
237c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
238c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
239c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
24036ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
241c58ca2e3SHong Zhang   PetscInt       nextb;
242c58ca2e3SHong Zhang 
243c58ca2e3SHong Zhang   PetscFunctionBegin;
244c58ca2e3SHong Zhang   /* clean old values in C */
245c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
246c58ca2e3SHong Zhang   /* Traverse A row-wise. */
247c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
248c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
249519eb980SHong Zhang   for (i=0; i<am; i++) {
250519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
251519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
252519eb980SHong Zhang     for (j=0; j<anzi; j++) {
253519eb980SHong Zhang       brow = aj[j];
254519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
255519eb980SHong Zhang       bjj  = bj + bi[brow];
256519eb980SHong Zhang       baj  = ba + bi[brow];
257519eb980SHong Zhang       /* perform sparse axpy */
25836ec6d2dSHong Zhang       valtmp = aa[j];
259c124e916SHong Zhang       nextb  = 0;
260c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
261c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
26236ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
263c124e916SHong Zhang         }
264d50806bdSBarry Smith       }
265d50806bdSBarry Smith       flops += 2*bnzi;
266d50806bdSBarry Smith     }
267519eb980SHong Zhang     aj += anzi; aa += anzi;
268519eb980SHong Zhang     cj += cnzi; ca += cnzi;
269d50806bdSBarry Smith   }
270c58ca2e3SHong Zhang 
271716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
272716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
273d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
274d50806bdSBarry Smith   PetscFunctionReturn(0);
275d50806bdSBarry Smith }
276bc011b1eSHong Zhang 
277e9e4536cSHong Zhang #undef __FUNCT__
2783c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2793c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
28025296bd5SBarry Smith {
28125296bd5SBarry Smith   PetscErrorCode     ierr;
28225296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
28325296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2843c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
28525296bd5SBarry Smith   MatScalar          *ca;
28625296bd5SBarry Smith   PetscReal          afill;
287eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
288eca6b86aSHong Zhang   PetscTable         ta;
2890298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
29025296bd5SBarry Smith 
29125296bd5SBarry Smith   PetscFunctionBegin;
2923c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2933c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2943c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
295854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
2963c50cad2SHong Zhang   ci[0] = 0;
2973c50cad2SHong Zhang 
2983c50cad2SHong Zhang   /* create and initialize a linked list */
299*c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
300*c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
301eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
302eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
303eca6b86aSHong Zhang 
304eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr);
3053c50cad2SHong Zhang 
3063c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
307f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
3083c50cad2SHong Zhang   current_space = free_space;
3093c50cad2SHong Zhang 
3103c50cad2SHong Zhang   /* Determine ci and cj */
3113c50cad2SHong Zhang   for (i=0; i<am; i++) {
3123c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3133c50cad2SHong Zhang     aj   = a->j + ai[i];
3143c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3153c50cad2SHong Zhang       brow = aj[j];
3163c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3173c50cad2SHong Zhang       bj   = b->j + bi[brow];
3183c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3193c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3203c50cad2SHong Zhang     }
3213c50cad2SHong Zhang     cnzi = lnk[1];
3223c50cad2SHong Zhang 
3233c50cad2SHong Zhang     /* If free space is not available, make more free space */
3243c50cad2SHong Zhang     /* Double the amount of total space in the list */
3253c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
326f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
3273c50cad2SHong Zhang       ndouble++;
3283c50cad2SHong Zhang     }
3293c50cad2SHong Zhang 
3303c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3313c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3322205254eSKarl Rupp 
3333c50cad2SHong Zhang     current_space->array           += cnzi;
3343c50cad2SHong Zhang     current_space->local_used      += cnzi;
3353c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3362205254eSKarl Rupp 
3373c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3383c50cad2SHong Zhang   }
3393c50cad2SHong Zhang 
3403c50cad2SHong Zhang   /* Column indices are in the list of free space */
3413c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3423c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
343854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
3443c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3453c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
34625296bd5SBarry Smith 
34725296bd5SBarry Smith   /* Allocate space for ca */
348854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
34925296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
35025296bd5SBarry Smith 
35125296bd5SBarry Smith   /* put together the new symbolic matrix */
352ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
35333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
35425296bd5SBarry Smith 
35525296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
35625296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
35725296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
35825296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
35925296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
36025296bd5SBarry Smith   c->nonew   = 0;
3612205254eSKarl Rupp 
36225296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
36325296bd5SBarry Smith 
36425296bd5SBarry Smith   /* set MatInfo */
36525296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
36625296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
36725296bd5SBarry Smith   c->maxnz                     = ci[am];
36825296bd5SBarry Smith   c->nz                        = ci[am];
3693c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
37025296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
37125296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
37225296bd5SBarry Smith 
37325296bd5SBarry Smith #if defined(PETSC_USE_INFO)
37425296bd5SBarry Smith   if (ci[am]) {
37557622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
37657622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
37725296bd5SBarry Smith   } else {
37825296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
37925296bd5SBarry Smith   }
38025296bd5SBarry Smith #endif
38125296bd5SBarry Smith   PetscFunctionReturn(0);
38225296bd5SBarry Smith }
38325296bd5SBarry Smith 
38425296bd5SBarry Smith 
38525296bd5SBarry Smith #undef __FUNCT__
38625023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
38725023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
388e9e4536cSHong Zhang {
389e9e4536cSHong Zhang   PetscErrorCode     ierr;
390e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
391bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
39225c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
393e9e4536cSHong Zhang   MatScalar          *ca;
394e9e4536cSHong Zhang   PetscReal          afill;
395eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
396eca6b86aSHong Zhang   PetscTable         ta;
3970298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
398e9e4536cSHong Zhang 
399e9e4536cSHong Zhang   PetscFunctionBegin;
400bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
401bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
402bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
403854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
404bd958071SHong Zhang   ci[0] = 0;
405bd958071SHong Zhang 
406bd958071SHong Zhang   /* create and initialize a linked list */
407*c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
408*c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
409eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
410eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
411eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
412bd958071SHong Zhang 
413bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
414f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
415bd958071SHong Zhang   current_space = free_space;
416bd958071SHong Zhang 
417bd958071SHong Zhang   /* Determine ci and cj */
418bd958071SHong Zhang   for (i=0; i<am; i++) {
419bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
420bd958071SHong Zhang     aj   = a->j + ai[i];
421bd958071SHong Zhang     for (j=0; j<anzi; j++) {
422bd958071SHong Zhang       brow = aj[j];
423bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
424bd958071SHong Zhang       bj   = b->j + bi[brow];
425bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
426bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
427bd958071SHong Zhang     }
428bd958071SHong Zhang     cnzi = lnk[0];
429bd958071SHong Zhang 
430bd958071SHong Zhang     /* If free space is not available, make more free space */
431bd958071SHong Zhang     /* Double the amount of total space in the list */
432bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
433f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
434bd958071SHong Zhang       ndouble++;
435bd958071SHong Zhang     }
436bd958071SHong Zhang 
437bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
438bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4392205254eSKarl Rupp 
440bd958071SHong Zhang     current_space->array           += cnzi;
441bd958071SHong Zhang     current_space->local_used      += cnzi;
442bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4432205254eSKarl Rupp 
444bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
445bd958071SHong Zhang   }
446bd958071SHong Zhang 
447bd958071SHong Zhang   /* Column indices are in the list of free space */
448bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
449bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
450854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
451bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
452bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
453e9e4536cSHong Zhang 
454e9e4536cSHong Zhang   /* Allocate space for ca */
455bd958071SHong Zhang   /*-----------------------*/
456854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
457e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
458e9e4536cSHong Zhang 
459e9e4536cSHong Zhang   /* put together the new symbolic matrix */
460ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
46133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
462e9e4536cSHong Zhang 
463e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
464e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
465e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
466e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
467e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
468e9e4536cSHong Zhang   c->nonew   = 0;
4692205254eSKarl Rupp 
47025023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
471e9e4536cSHong Zhang 
472e9e4536cSHong Zhang   /* set MatInfo */
473e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
474e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
475e9e4536cSHong Zhang   c->maxnz                     = ci[am];
476e9e4536cSHong Zhang   c->nz                        = ci[am];
477bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
478e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
479e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
480e9e4536cSHong Zhang 
481e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
482e9e4536cSHong Zhang   if (ci[am]) {
48357622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
48457622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
485e9e4536cSHong Zhang   } else {
486e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
487e9e4536cSHong Zhang   }
488e9e4536cSHong Zhang #endif
489e9e4536cSHong Zhang   PetscFunctionReturn(0);
490e9e4536cSHong Zhang }
491e9e4536cSHong Zhang 
4920ced3a2bSJed Brown #undef __FUNCT__
4930ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
4940ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
4950ced3a2bSJed Brown {
4960ced3a2bSJed Brown   PetscErrorCode     ierr;
4970ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4980ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
4990ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
5000ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
5010ced3a2bSJed Brown   PetscReal          afill;
5020ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
5030298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
5040ced3a2bSJed Brown   PetscHeap          h;
5050ced3a2bSJed Brown 
5060ced3a2bSJed Brown   PetscFunctionBegin;
507cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
5080ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
5090ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
510854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
5110ced3a2bSJed Brown   ci[0] = 0;
5120ced3a2bSJed Brown 
5130ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
514f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
5150ced3a2bSJed Brown   current_space = free_space;
5160ced3a2bSJed Brown 
5170ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
518785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
5190ced3a2bSJed Brown 
5200ced3a2bSJed Brown   /* Determine ci and cj */
5210ced3a2bSJed Brown   for (i=0; i<am; i++) {
5220ced3a2bSJed Brown     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
5230ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5240ced3a2bSJed Brown     ci[i+1] = ci[i];
5250ced3a2bSJed Brown     /* Populate the min heap */
5260ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5270ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5280ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5290ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5300ced3a2bSJed Brown       }
5310ced3a2bSJed Brown     }
5320ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5330ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5340ced3a2bSJed Brown     while (j >= 0) {
5350ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
536f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
5370ced3a2bSJed Brown         ndouble++;
5380ced3a2bSJed Brown       }
5390ced3a2bSJed Brown       *(current_space->array++) = col;
5400ced3a2bSJed Brown       current_space->local_used++;
5410ced3a2bSJed Brown       current_space->local_remaining--;
5420ced3a2bSJed Brown       ci[i+1]++;
5430ced3a2bSJed Brown 
5440ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5450ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5460ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5470ced3a2bSJed Brown         PetscInt j2,col2;
5480ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5490ced3a2bSJed Brown         if (col2 != col) break;
5500ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5510ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5520ced3a2bSJed Brown       }
5530ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5540ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5550ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5560ced3a2bSJed Brown     }
5570ced3a2bSJed Brown   }
5580ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5590ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5600ced3a2bSJed Brown 
5610ced3a2bSJed Brown   /* Column indices are in the list of free space */
5620ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5630ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
564785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
5650ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5660ced3a2bSJed Brown 
5670ced3a2bSJed Brown   /* put together the new symbolic matrix */
568ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
56933d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
5700ced3a2bSJed Brown 
5710ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5720ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5730ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5740ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5750ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5760ced3a2bSJed Brown   c->nonew   = 0;
57726fbe8dcSKarl Rupp 
57889d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5790ced3a2bSJed Brown 
5800ced3a2bSJed Brown   /* set MatInfo */
5810ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5820ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5830ced3a2bSJed Brown   c->maxnz                     = ci[am];
5840ced3a2bSJed Brown   c->nz                        = ci[am];
5850ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5860ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5870ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5880ced3a2bSJed Brown 
5890ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5900ced3a2bSJed Brown   if (ci[am]) {
59157622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
59257622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
5930ced3a2bSJed Brown   } else {
5940ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
5950ced3a2bSJed Brown   }
5960ced3a2bSJed Brown #endif
5970ced3a2bSJed Brown   PetscFunctionReturn(0);
5980ced3a2bSJed Brown }
599e9e4536cSHong Zhang 
6008a07c6f1SJed Brown #undef __FUNCT__
6018a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
6028a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
6038a07c6f1SJed Brown {
6048a07c6f1SJed Brown   PetscErrorCode     ierr;
6058a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6068a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
6078a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6088a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6098a07c6f1SJed Brown   PetscReal          afill;
6108a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6110298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6128a07c6f1SJed Brown   PetscHeap          h;
6138a07c6f1SJed Brown   PetscBT            bt;
6148a07c6f1SJed Brown 
6158a07c6f1SJed Brown   PetscFunctionBegin;
616cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6178a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6188a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
619854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
6208a07c6f1SJed Brown   ci[0] = 0;
6218a07c6f1SJed Brown 
6228a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
623f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
6242205254eSKarl Rupp 
6258a07c6f1SJed Brown   current_space = free_space;
6268a07c6f1SJed Brown 
6278a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
628785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
6298a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6308a07c6f1SJed Brown 
6318a07c6f1SJed Brown   /* Determine ci and cj */
6328a07c6f1SJed Brown   for (i=0; i<am; i++) {
6338a07c6f1SJed Brown     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
6348a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6358a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6368a07c6f1SJed Brown     ci[i+1] = ci[i];
6378a07c6f1SJed Brown     /* Populate the min heap */
6388a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6398a07c6f1SJed Brown       PetscInt brow = acol[j];
6408a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6418a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6428a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6438a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6448a07c6f1SJed Brown           bb[j]++;
6458a07c6f1SJed Brown           break;
6468a07c6f1SJed Brown         }
6478a07c6f1SJed Brown       }
6488a07c6f1SJed Brown     }
6498a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6508a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6518a07c6f1SJed Brown     while (j >= 0) {
6528a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6530298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
654f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
6558a07c6f1SJed Brown         ndouble++;
6568a07c6f1SJed Brown       }
6578a07c6f1SJed Brown       *(current_space->array++) = col;
6588a07c6f1SJed Brown       current_space->local_used++;
6598a07c6f1SJed Brown       current_space->local_remaining--;
6608a07c6f1SJed Brown       ci[i+1]++;
6618a07c6f1SJed Brown 
6628a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6638a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6648a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6658a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6668a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6678a07c6f1SJed Brown           bb[j]++;
6688a07c6f1SJed Brown           break;
6698a07c6f1SJed Brown         }
6708a07c6f1SJed Brown       }
6718a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6728a07c6f1SJed Brown     }
6738a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6748a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6758a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6768a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6778a07c6f1SJed Brown     }
6788a07c6f1SJed Brown   }
6798a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6808a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6818a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6828a07c6f1SJed Brown 
6838a07c6f1SJed Brown   /* Column indices are in the list of free space */
6848a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6858a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
686785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
6878a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6888a07c6f1SJed Brown 
6898a07c6f1SJed Brown   /* put together the new symbolic matrix */
690ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
69133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
6928a07c6f1SJed Brown 
6938a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6948a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6958a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
6968a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
6978a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
6988a07c6f1SJed Brown   c->nonew   = 0;
69926fbe8dcSKarl Rupp 
70089d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
7018a07c6f1SJed Brown 
7028a07c6f1SJed Brown   /* set MatInfo */
7038a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
7048a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
7058a07c6f1SJed Brown   c->maxnz                     = ci[am];
7068a07c6f1SJed Brown   c->nz                        = ci[am];
7078a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7088a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7098a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7108a07c6f1SJed Brown 
7118a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7128a07c6f1SJed Brown   if (ci[am]) {
71357622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
71457622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
7158a07c6f1SJed Brown   } else {
7168a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7178a07c6f1SJed Brown   }
7188a07c6f1SJed Brown #endif
7198a07c6f1SJed Brown   PetscFunctionReturn(0);
7208a07c6f1SJed Brown }
7218a07c6f1SJed Brown 
722cd093f1dSJed Brown #undef __FUNCT__
72358cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
724cd093f1dSJed Brown /* concatenate unique entries and then sort */
72558cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
726cd093f1dSJed Brown {
727cd093f1dSJed Brown   PetscErrorCode     ierr;
728cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
729cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
730cd093f1dSJed Brown   PetscInt           *ci,*cj;
731cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
732cd093f1dSJed Brown   PetscReal          afill;
733cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
734cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
735cd093f1dSJed Brown   char               *seen;
736cd093f1dSJed Brown 
737cd093f1dSJed Brown   PetscFunctionBegin;
738854ce69bSBarry Smith   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
739cd093f1dSJed Brown   ci[0] = 0;
740cd093f1dSJed Brown 
741cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
742cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
743cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
744785e854fSJed Brown   ierr = PetscMalloc1(bn,&seen);CHKERRQ(ierr);
745cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
746cd093f1dSJed Brown 
747cd093f1dSJed Brown   /* Determine ci and cj */
748cd093f1dSJed Brown   for (i=0; i<am; i++) {
749cd093f1dSJed Brown     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
750cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
751cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
752cd093f1dSJed Brown     /* Pack segrow */
753cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
754cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
755cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
756cd093f1dSJed Brown         PetscInt bcol = bj[k];
757cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
758cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
759cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
760cd093f1dSJed Brown           *slot = bcol;
761cd093f1dSJed Brown           seen[bcol] = 1;
762cd093f1dSJed Brown           packlen++;
763cd093f1dSJed Brown         }
764cd093f1dSJed Brown       }
765cd093f1dSJed Brown     }
766cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
767cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
768cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
769cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
770cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
771cd093f1dSJed Brown   }
772cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
773cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
774cd093f1dSJed Brown 
775cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
776cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
777cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
778cd093f1dSJed Brown 
779cd093f1dSJed Brown   /* put together the new symbolic matrix */
780cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
78133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
782cd093f1dSJed Brown 
783cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
784cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
785cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
786cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
787cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
788cd093f1dSJed Brown   c->nonew   = 0;
789cd093f1dSJed Brown 
790cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
791cd093f1dSJed Brown 
792cd093f1dSJed Brown   /* set MatInfo */
793cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
794cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
795cd093f1dSJed Brown   c->maxnz                     = ci[am];
796cd093f1dSJed Brown   c->nz                        = ci[am];
797cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
798cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
799cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
800cd093f1dSJed Brown 
801cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
802cd093f1dSJed Brown   if (ci[am]) {
80357622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
80457622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
805cd093f1dSJed Brown   } else {
806cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
807cd093f1dSJed Brown   }
808cd093f1dSJed Brown #endif
809cd093f1dSJed Brown   PetscFunctionReturn(0);
810cd093f1dSJed Brown }
811cd093f1dSJed Brown 
812d2853540SHong Zhang /* This routine is not used. Should be removed! */
813bc011b1eSHong Zhang #undef __FUNCT__
8146fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
8156fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8165df89d91SHong Zhang {
817bc011b1eSHong Zhang   PetscErrorCode ierr;
818bc011b1eSHong Zhang 
819bc011b1eSHong Zhang   PetscFunctionBegin;
820bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
8213ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
8226fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8233ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
824bc011b1eSHong Zhang   }
8253ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
8266fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8273ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
828bc011b1eSHong Zhang   PetscFunctionReturn(0);
829bc011b1eSHong Zhang }
830bc011b1eSHong Zhang 
831bc011b1eSHong Zhang #undef __FUNCT__
8322128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
8332128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8342128a86cSHong Zhang {
8352128a86cSHong Zhang   PetscErrorCode      ierr;
8364c7df5ccSHong Zhang   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
83740192850SHong Zhang   Mat_MatMatTransMult *abt=a->abt;
8382128a86cSHong Zhang 
8392128a86cSHong Zhang   PetscFunctionBegin;
84040192850SHong Zhang   ierr = (abt->destroy)(A);CHKERRQ(ierr);
84140192850SHong Zhang   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
84240192850SHong Zhang   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
84340192850SHong Zhang   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
84440192850SHong Zhang   ierr = PetscFree(abt);CHKERRQ(ierr);
8452128a86cSHong Zhang   PetscFunctionReturn(0);
8462128a86cSHong Zhang }
8472128a86cSHong Zhang 
8482128a86cSHong Zhang #undef __FUNCT__
8496fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
8506fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
851bc011b1eSHong Zhang {
8525df89d91SHong Zhang   PetscErrorCode      ierr;
85381d82fe4SHong Zhang   Mat                 Bt;
85481d82fe4SHong Zhang   PetscInt            *bti,*btj;
85540192850SHong Zhang   Mat_MatMatTransMult *abt;
8564c7df5ccSHong Zhang   Mat_SeqAIJ          *c;
857d2853540SHong Zhang 
85881d82fe4SHong Zhang   PetscFunctionBegin;
85981d82fe4SHong Zhang   /* create symbolic Bt */
86081d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8610298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
86233d57670SJed Brown   ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
86381d82fe4SHong Zhang 
86481d82fe4SHong Zhang   /* get symbolic C=A*Bt */
86581d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
86681d82fe4SHong Zhang 
8672128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
868b00a9115SJed Brown   ierr   = PetscNew(&abt);CHKERRQ(ierr);
8694c7df5ccSHong Zhang   c      = (Mat_SeqAIJ*)(*C)->data;
87040192850SHong Zhang   c->abt = abt;
8712128a86cSHong Zhang 
87240192850SHong Zhang   abt->usecoloring = PETSC_FALSE;
87340192850SHong Zhang   abt->destroy     = (*C)->ops->destroy;
8742128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8752128a86cSHong Zhang 
876c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr);
87740192850SHong Zhang   if (abt->usecoloring) {
878b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
879b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
880335efc43SPeter Brune     MatColoring          coloring;
8812128a86cSHong Zhang     ISColoring           iscoloring;
8822128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
8834d478ae7SHong Zhang     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
8844d478ae7SHong Zhang     /* inode causes memory problem, don't know why */
8854d478ae7SHong Zhang     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
886e8354b3eSHong Zhang 
887335efc43SPeter Brune     ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr);
888335efc43SPeter Brune     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
889335efc43SPeter Brune     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
890335efc43SPeter Brune     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
891335efc43SPeter Brune     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
892335efc43SPeter Brune     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
893b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
8942205254eSKarl Rupp 
89540192850SHong Zhang     abt->matcoloring = matcoloring;
8962205254eSKarl Rupp 
8972128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
8982128a86cSHong Zhang 
8992128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
9002128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
9012128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9022128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
9030298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9042205254eSKarl Rupp 
9052128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
90640192850SHong Zhang     abt->Bt_den   = Bt_dense;
9072128a86cSHong Zhang 
9082128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9092128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9102128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9110298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9122205254eSKarl Rupp 
9132128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
91440192850SHong Zhang     abt->ABt_den  = C_dense;
915f94ccd6cSHong Zhang 
916f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9171ce71dffSSatish Balay     {
918f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
919c40ebe3bSHong Zhang       ierr = PetscInfo7(*C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr);
9201ce71dffSSatish Balay     }
921f94ccd6cSHong Zhang #endif
9222128a86cSHong Zhang   }
923e99be685SHong Zhang   /* clean up */
924e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
925e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9265df89d91SHong Zhang   PetscFunctionReturn(0);
9275df89d91SHong Zhang }
9285df89d91SHong Zhang 
9295df89d91SHong Zhang #undef __FUNCT__
9306fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9316fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9325df89d91SHong Zhang {
9335df89d91SHong Zhang   PetscErrorCode      ierr;
9345df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
935e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
93681d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9375df89d91SHong Zhang   PetscLogDouble      flops=0.0;
938aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
93940192850SHong Zhang   Mat_MatMatTransMult *abt = c->abt;
9405df89d91SHong Zhang 
9415df89d91SHong Zhang   PetscFunctionBegin;
94258ed3ceaSHong Zhang   /* clear old values in C */
94358ed3ceaSHong Zhang   if (!c->a) {
944854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
94558ed3ceaSHong Zhang     c->a      = ca;
94658ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
94758ed3ceaSHong Zhang   } else {
94858ed3ceaSHong Zhang     ca =  c->a;
94958ed3ceaSHong Zhang   }
95058ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
95158ed3ceaSHong Zhang 
95240192850SHong Zhang   if (abt->usecoloring) {
95340192850SHong Zhang     MatTransposeColoring matcoloring = abt->matcoloring;
95440192850SHong Zhang     Mat                  Bt_dense,C_dense = abt->ABt_den;
955c8db22f6SHong Zhang 
956b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
95740192850SHong Zhang     Bt_dense = abt->Bt_den;
958b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
959c8db22f6SHong Zhang 
960c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9612128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
962c8db22f6SHong Zhang 
9632128a86cSHong Zhang     /* Recover C from C_dense */
964b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
965c8db22f6SHong Zhang     PetscFunctionReturn(0);
966c8db22f6SHong Zhang   }
967c8db22f6SHong Zhang 
9681fa1209cSHong Zhang   for (i=0; i<cm; i++) {
96981d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
97081d82fe4SHong Zhang     acol = aj + ai[i];
9716973769bSHong Zhang     aval = aa + ai[i];
9721fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9731fa1209cSHong Zhang     ccol = cj + ci[i];
9746973769bSHong Zhang     cval = ca + ci[i];
9751fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
97681d82fe4SHong Zhang       brow = ccol[j];
97781d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
97881d82fe4SHong Zhang       bcol = bj + bi[brow];
9796973769bSHong Zhang       bval = ba + bi[brow];
9806973769bSHong Zhang 
98181d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
98281d82fe4SHong Zhang       nexta = 0; nextb = 0;
98381d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
9847b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
98581d82fe4SHong Zhang         if (nexta == anzi) break;
9867b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
98781d82fe4SHong Zhang         if (nextb == bnzj) break;
98881d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
9896973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
99081d82fe4SHong Zhang           nexta++; nextb++;
99181d82fe4SHong Zhang           flops += 2;
9921fa1209cSHong Zhang         }
9931fa1209cSHong Zhang       }
99481d82fe4SHong Zhang     }
99581d82fe4SHong Zhang   }
99681d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99781d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99881d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
9995df89d91SHong Zhang   PetscFunctionReturn(0);
10005df89d91SHong Zhang }
10015df89d91SHong Zhang 
10025df89d91SHong Zhang #undef __FUNCT__
100375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
10040adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10050adebc6cSBarry Smith {
10065df89d91SHong Zhang   PetscErrorCode ierr;
10075df89d91SHong Zhang 
10085df89d91SHong Zhang   PetscFunctionBegin;
10095df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
101007706670SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
101175648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
101207706670SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
10135df89d91SHong Zhang   }
101407706670SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
101575648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
101607706670SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
10175df89d91SHong Zhang   PetscFunctionReturn(0);
10185df89d91SHong Zhang }
10195df89d91SHong Zhang 
10205df89d91SHong Zhang #undef __FUNCT__
102175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
102275648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10235df89d91SHong Zhang {
1024bc011b1eSHong Zhang   PetscErrorCode ierr;
1025bc011b1eSHong Zhang   Mat            At;
102638baddfdSBarry Smith   PetscInt       *ati,*atj;
1027bc011b1eSHong Zhang 
1028bc011b1eSHong Zhang   PetscFunctionBegin;
1029bc011b1eSHong Zhang   /* create symbolic At */
1030bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10310298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
103233d57670SJed Brown   ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1033bc011b1eSHong Zhang 
1034bc011b1eSHong Zhang   /* get symbolic C=At*B */
1035bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1036bc011b1eSHong Zhang 
1037bc011b1eSHong Zhang   /* clean up */
10386bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1039bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1040bc011b1eSHong Zhang   PetscFunctionReturn(0);
1041bc011b1eSHong Zhang }
1042bc011b1eSHong Zhang 
1043bc011b1eSHong Zhang #undef __FUNCT__
104475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
104575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1046bc011b1eSHong Zhang {
1047bc011b1eSHong Zhang   PetscErrorCode ierr;
10480fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1049d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1050d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1051d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1052aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1053bc011b1eSHong Zhang 
1054bc011b1eSHong Zhang   PetscFunctionBegin;
1055aa1aec99SHong Zhang   if (!c->a) {
1056854ce69bSBarry Smith     ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
10572205254eSKarl Rupp 
1058aa1aec99SHong Zhang     c->a      = ca;
1059aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1060aa1aec99SHong Zhang   } else {
1061aa1aec99SHong Zhang     ca = c->a;
1062aa1aec99SHong Zhang   }
1063bc011b1eSHong Zhang   /* clear old values in C */
1064bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1065bc011b1eSHong Zhang 
1066bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1067bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1068bc011b1eSHong Zhang     bj   = b->j + bi[i];
1069bc011b1eSHong Zhang     ba   = b->a + bi[i];
1070bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1071bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1072bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1073bc011b1eSHong Zhang       nextb = 0;
10740fbc74f4SHong Zhang       crow  = *aj++;
1075bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1076bc011b1eSHong Zhang       caj   = ca + ci[crow];
1077bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1078bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
10790fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
10800fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1081bc011b1eSHong Zhang           nextb++;
1082bc011b1eSHong Zhang         }
1083bc011b1eSHong Zhang       }
1084bc011b1eSHong Zhang       flops += 2*bnzi;
10850fbc74f4SHong Zhang       aa++;
1086bc011b1eSHong Zhang     }
1087bc011b1eSHong Zhang   }
1088bc011b1eSHong Zhang 
1089bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1090bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1091bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1092bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1093bc011b1eSHong Zhang   PetscFunctionReturn(0);
1094bc011b1eSHong Zhang }
10959123193aSHong Zhang 
10969123193aSHong Zhang #undef __FUNCT__
10979123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
1098150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10999123193aSHong Zhang {
11009123193aSHong Zhang   PetscErrorCode ierr;
11019123193aSHong Zhang 
11029123193aSHong Zhang   PetscFunctionBegin;
11039123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
11043ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11059123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11063ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11079123193aSHong Zhang   }
11083ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11099123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11104614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11119123193aSHong Zhang   PetscFunctionReturn(0);
11129123193aSHong Zhang }
1113edd81eecSMatthew Knepley 
11149123193aSHong Zhang #undef __FUNCT__
11159123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11169123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11179123193aSHong Zhang {
11189123193aSHong Zhang   PetscErrorCode ierr;
11199123193aSHong Zhang 
11209123193aSHong Zhang   PetscFunctionBegin;
11215a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11222205254eSKarl Rupp 
1123d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11249123193aSHong Zhang   PetscFunctionReturn(0);
11259123193aSHong Zhang }
11269123193aSHong Zhang 
11279123193aSHong Zhang #undef __FUNCT__
11289123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11299123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11309123193aSHong Zhang {
1131f32d5d43SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1132612438f1SStefano Zampini   Mat_SeqDense      *bd = (Mat_SeqDense*)B->data;
11339123193aSHong Zhang   PetscErrorCode    ierr;
1134daf57211SHong Zhang   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1135daf57211SHong Zhang   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1136daf57211SHong Zhang   const PetscInt    *aj;
1137612438f1SStefano Zampini   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1138daf57211SHong Zhang   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
11399123193aSHong Zhang 
11409123193aSHong Zhang   PetscFunctionBegin;
1141f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1142612438f1SStefano Zampini   if (B->rmap->n != 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,B->rmap->n);
1143e32f2f54SBarry 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);
1144e32f2f54SBarry 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);
1145612438f1SStefano Zampini   b = bd->v;
11468c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1147f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1148730858b9SHong Zhang   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1149f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1150f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1151f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1152f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1153f32d5d43SBarry Smith       aj = a->j + a->i[i];
1154f32d5d43SBarry Smith       aa = a->a + a->i[i];
1155f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1156730858b9SHong Zhang         aatmp = aa[j]; ajtmp = aj[j];
1157730858b9SHong Zhang         r1 += aatmp*b1[ajtmp];
1158730858b9SHong Zhang         r2 += aatmp*b2[ajtmp];
1159730858b9SHong Zhang         r3 += aatmp*b3[ajtmp];
1160730858b9SHong Zhang         r4 += aatmp*b4[ajtmp];
11619123193aSHong Zhang       }
1162730858b9SHong Zhang       c1[i] = r1;
1163730858b9SHong Zhang       c2[i] = r2;
1164730858b9SHong Zhang       c3[i] = r3;
1165730858b9SHong Zhang       c4[i] = r4;
1166f32d5d43SBarry Smith     }
1167730858b9SHong Zhang     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1168730858b9SHong Zhang     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1169f32d5d43SBarry Smith   }
1170f32d5d43SBarry Smith   for (; col<cn; col++) {   /* over extra columns of C */
1171f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1172f32d5d43SBarry Smith       r1 = 0.0;
1173f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1174f32d5d43SBarry Smith       aj = a->j + a->i[i];
1175f32d5d43SBarry Smith       aa = a->a + a->i[i];
1176f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1177730858b9SHong Zhang         r1 += aa[j]*b1[aj[j]];
1178f32d5d43SBarry Smith       }
1179730858b9SHong Zhang       c1[i] = r1;
1180f32d5d43SBarry Smith     }
1181f32d5d43SBarry Smith     b1 += bm;
1182730858b9SHong Zhang     c1 += am;
1183f32d5d43SBarry Smith   }
1184dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
11858c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1186f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1187f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1188f32d5d43SBarry Smith   PetscFunctionReturn(0);
1189f32d5d43SBarry Smith }
1190f32d5d43SBarry Smith 
1191f32d5d43SBarry Smith /*
11924324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1193f32d5d43SBarry Smith */
1194f32d5d43SBarry Smith #undef __FUNCT__
1195f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1196f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1197f32d5d43SBarry Smith {
1198f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
119990f5ea3eSStefano Zampini   Mat_SeqDense   *bd = (Mat_SeqDense*)B->data;
1200f32d5d43SBarry Smith   PetscErrorCode ierr;
1201dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1202dd6ea824SBarry Smith   MatScalar      *aa;
120390f5ea3eSStefano Zampini   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
12044324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1205f32d5d43SBarry Smith 
1206f32d5d43SBarry Smith   PetscFunctionBegin;
1207f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
120890f5ea3eSStefano Zampini   b = bd->v;
12098c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1210f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12114324174eSBarry Smith 
12124324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12134324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12144324174eSBarry Smith       colam = col*am;
12154324174eSBarry Smith       arm   = a->compressedrow.nrows;
12164324174eSBarry Smith       ii    = a->compressedrow.i;
12174324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12184324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12194324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12204324174eSBarry Smith         n  = ii[i+1] - ii[i];
12214324174eSBarry Smith         aj = a->j + ii[i];
12224324174eSBarry Smith         aa = a->a + ii[i];
12234324174eSBarry Smith         for (j=0; j<n; j++) {
12244324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12254324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12264324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12274324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12284324174eSBarry Smith         }
12294324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12304324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12314324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12324324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12334324174eSBarry Smith       }
12344324174eSBarry Smith       b1 += bm4;
12354324174eSBarry Smith       b2 += bm4;
12364324174eSBarry Smith       b3 += bm4;
12374324174eSBarry Smith       b4 += bm4;
12384324174eSBarry Smith     }
12394324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12404324174eSBarry Smith       colam = col*am;
12414324174eSBarry Smith       arm   = a->compressedrow.nrows;
12424324174eSBarry Smith       ii    = a->compressedrow.i;
12434324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12444324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12454324174eSBarry Smith         r1 = 0.0;
12464324174eSBarry Smith         n  = ii[i+1] - ii[i];
12474324174eSBarry Smith         aj = a->j + ii[i];
12484324174eSBarry Smith         aa = a->a + ii[i];
12494324174eSBarry Smith 
12504324174eSBarry Smith         for (j=0; j<n; j++) {
12514324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
12524324174eSBarry Smith         }
1253a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
12544324174eSBarry Smith       }
12554324174eSBarry Smith       b1 += bm;
12564324174eSBarry Smith     }
12574324174eSBarry Smith   } else {
1258f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1259f32d5d43SBarry Smith       colam = col*am;
1260f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1261f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1262f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1263f32d5d43SBarry Smith         aj = a->j + a->i[i];
1264f32d5d43SBarry Smith         aa = a->a + a->i[i];
1265f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1266f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1267f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1268f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1269f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1270f32d5d43SBarry Smith         }
1271f32d5d43SBarry Smith         c[colam + i]       += r1;
1272f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1273f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1274f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1275f32d5d43SBarry Smith       }
1276f32d5d43SBarry Smith       b1 += bm4;
1277f32d5d43SBarry Smith       b2 += bm4;
1278f32d5d43SBarry Smith       b3 += bm4;
1279f32d5d43SBarry Smith       b4 += bm4;
1280f32d5d43SBarry Smith     }
1281f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1282a2ea699eSBarry Smith       colam = col*am;
1283f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1284f32d5d43SBarry Smith         r1 = 0.0;
1285f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1286f32d5d43SBarry Smith         aj = a->j + a->i[i];
1287f32d5d43SBarry Smith         aa = a->a + a->i[i];
1288f32d5d43SBarry Smith 
1289f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1290f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1291f32d5d43SBarry Smith         }
1292a2ea699eSBarry Smith         c[colam + i] += r1;
1293f32d5d43SBarry Smith       }
1294f32d5d43SBarry Smith       b1 += bm;
1295f32d5d43SBarry Smith     }
12964324174eSBarry Smith   }
1297dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
12988c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
12999123193aSHong Zhang   PetscFunctionReturn(0);
13009123193aSHong Zhang }
1301b1683b59SHong Zhang 
1302b1683b59SHong Zhang #undef __FUNCT__
1303b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1304b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1305c8db22f6SHong Zhang {
1306c8db22f6SHong Zhang   PetscErrorCode ierr;
13072f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13082f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13092f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13102f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13112f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13122f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1313c8db22f6SHong Zhang 
1314c8db22f6SHong Zhang   PetscFunctionBegin;
13152f5f1f90SHong Zhang   btval_den=btdense->v;
13162f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13172f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13182f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13192f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1320d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13212f5f1f90SHong Zhang       btcol = bj + bi[col];
13222f5f1f90SHong Zhang       btval = ba + bi[col];
13232f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1324d2853540SHong Zhang       for (j=0; j<anz; j++) {
13252f5f1f90SHong Zhang         brow            = btcol[j];
13262f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1327c8db22f6SHong Zhang       }
1328c8db22f6SHong Zhang     }
13292f5f1f90SHong Zhang     btval_den += m;
1330c8db22f6SHong Zhang   }
1331c8db22f6SHong Zhang   PetscFunctionReturn(0);
1332c8db22f6SHong Zhang }
1333c8db22f6SHong Zhang 
1334c8db22f6SHong Zhang #undef __FUNCT__
1335b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1336b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13378972f759SHong Zhang {
13388972f759SHong Zhang   PetscErrorCode ierr;
1339b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
1340077f23c2SHong Zhang   PetscScalar    *ca_den,*ca_den_ptr,*ca=csp->a;
1341f99a636bSHong Zhang   PetscInt       k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1342e88777f2SHong Zhang   PetscInt       brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1343077f23c2SHong Zhang   PetscInt       nrows,*row,*idx;
1344077f23c2SHong Zhang   PetscInt       *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
13458972f759SHong Zhang 
13468972f759SHong Zhang   PetscFunctionBegin;
1347a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1348f99a636bSHong Zhang 
1349077f23c2SHong Zhang   if (brows > 0) {
1350077f23c2SHong Zhang     PetscInt *lstart,row_end,row_start;
1351980ae229SHong Zhang     lstart = matcoloring->lstart;
1352eeb4fd02SHong Zhang     ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr);
1353eeb4fd02SHong Zhang 
1354077f23c2SHong Zhang     row_end = brows;
1355eeb4fd02SHong Zhang     if (row_end > m) row_end = m;
1356077f23c2SHong Zhang     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1357077f23c2SHong Zhang       ca_den_ptr = ca_den;
1358980ae229SHong Zhang       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1359eeb4fd02SHong Zhang         nrows = matcoloring->nrows[k];
1360eeb4fd02SHong Zhang         row   = rows  + colorforrow[k];
1361eeb4fd02SHong Zhang         idx   = den2sp + colorforrow[k];
1362eeb4fd02SHong Zhang         for (l=lstart[k]; l<nrows; l++) {
1363eeb4fd02SHong Zhang           if (row[l] >= row_end) {
1364eeb4fd02SHong Zhang             lstart[k] = l;
1365eeb4fd02SHong Zhang             break;
1366eeb4fd02SHong Zhang           } else {
1367077f23c2SHong Zhang             ca[idx[l]] = ca_den_ptr[row[l]];
1368eeb4fd02SHong Zhang           }
1369eeb4fd02SHong Zhang         }
1370077f23c2SHong Zhang         ca_den_ptr += m;
1371eeb4fd02SHong Zhang       }
1372077f23c2SHong Zhang       row_end += brows;
1373eeb4fd02SHong Zhang       if (row_end > m) row_end = m;
1374eeb4fd02SHong Zhang     }
1375077f23c2SHong Zhang   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1376077f23c2SHong Zhang     ca_den_ptr = ca_den;
1377077f23c2SHong Zhang     for (k=0; k<ncolors; k++) {
1378077f23c2SHong Zhang       nrows = matcoloring->nrows[k];
1379077f23c2SHong Zhang       row   = rows  + colorforrow[k];
1380077f23c2SHong Zhang       idx   = den2sp + colorforrow[k];
1381077f23c2SHong Zhang       for (l=0; l<nrows; l++) {
1382077f23c2SHong Zhang         ca[idx[l]] = ca_den_ptr[row[l]];
1383077f23c2SHong Zhang       }
1384077f23c2SHong Zhang       ca_den_ptr += m;
1385077f23c2SHong Zhang     }
1386f99a636bSHong Zhang   }
1387f99a636bSHong Zhang 
1388a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
1389f99a636bSHong Zhang #if defined(PETSC_USE_INFO)
1390077f23c2SHong Zhang   if (matcoloring->brows > 0) {
1391f99a636bSHong Zhang     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1392e88777f2SHong Zhang   } else {
1393077f23c2SHong Zhang     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1394e88777f2SHong Zhang   }
1395f99a636bSHong Zhang #endif
13968972f759SHong Zhang   PetscFunctionReturn(0);
13978972f759SHong Zhang }
13988972f759SHong Zhang 
13998972f759SHong Zhang #undef __FUNCT__
1400b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1401b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1402b1683b59SHong Zhang {
1403b1683b59SHong Zhang   PetscErrorCode ierr;
1404e88777f2SHong Zhang   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
14051a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1406b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1407b1683b59SHong Zhang   IS             *isa;
1408b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1409e88777f2SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1410e88777f2SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i,brows;
1411e88777f2SHong Zhang   PetscBool      flg;
1412b1683b59SHong Zhang 
1413b1683b59SHong Zhang   PetscFunctionBegin;
1414b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1415e99be685SHong Zhang 
14167ecccc15SHong Zhang   /* bs >1 is not being tested yet! */
1417e88777f2SHong Zhang   Nbs       = mat->cmap->N/bs;
1418b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1419e88777f2SHong Zhang   c->N      = Nbs;
1420e88777f2SHong Zhang   c->m      = c->M;
1421b1683b59SHong Zhang   c->rstart = 0;
1422e88777f2SHong Zhang   c->brows  = 100;
1423b1683b59SHong Zhang 
1424b1683b59SHong Zhang   c->ncolors = nis;
1425dcca6d9dSJed Brown   ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr);
1426854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr);
1427854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr);
1428e88777f2SHong Zhang 
1429e88777f2SHong Zhang   brows = c->brows;
1430c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1431e88777f2SHong Zhang   if (flg) c->brows = brows;
1432eeb4fd02SHong Zhang   if (brows > 0) {
1433854ce69bSBarry Smith     ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr);
1434e88777f2SHong Zhang   }
14352205254eSKarl Rupp 
1436d2853540SHong Zhang   colorforrow[0] = 0;
1437d2853540SHong Zhang   rows_i         = rows;
1438f99a636bSHong Zhang   den2sp_i       = den2sp;
1439b1683b59SHong Zhang 
1440854ce69bSBarry Smith   ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr);
1441854ce69bSBarry Smith   ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr);
14422205254eSKarl Rupp 
1443d2853540SHong Zhang   colorforcol[0] = 0;
1444d2853540SHong Zhang   columns_i      = columns;
1445d2853540SHong Zhang 
1446077f23c2SHong Zhang   /* get column-wise storage of mat */
1447077f23c2SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1448b1683b59SHong Zhang 
1449b28d80bdSHong Zhang   cm   = c->m;
1450854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr);
1451854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr);
1452f99a636bSHong Zhang   for (i=0; i<nis; i++) { /* loop over color */
1453b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1454b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
14552205254eSKarl Rupp 
1456b1683b59SHong Zhang     c->ncolumns[i] = n;
1457b1683b59SHong Zhang     if (n) {
1458d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1459b1683b59SHong Zhang     }
1460d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1461d2853540SHong Zhang     columns_i       += n;
1462b1683b59SHong Zhang 
1463b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1464e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1465e99be685SHong Zhang 
1466b7caf3d6SHong Zhang     for (j=0; j<n; j++) { /* loop over columns*/
1467b1683b59SHong Zhang       col     = is[j];
1468d2853540SHong Zhang       row_idx = cj + ci[col];
1469b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1470b7caf3d6SHong Zhang       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1471e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1472d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1473b1683b59SHong Zhang       }
1474b1683b59SHong Zhang     }
1475b1683b59SHong Zhang     /* count the number of hits */
1476b1683b59SHong Zhang     nrows = 0;
1477e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1478b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1479b1683b59SHong Zhang     }
1480b1683b59SHong Zhang     c->nrows[i]      = nrows;
1481d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1482d2853540SHong Zhang 
1483b1683b59SHong Zhang     nrows = 0;
1484b7caf3d6SHong Zhang     for (j=0; j<cm; j++) { /* loop over rows */
1485b1683b59SHong Zhang       if (rowhit[j]) {
1486d2853540SHong Zhang         rows_i[nrows]   = j;
148712b89a43SHong Zhang         den2sp_i[nrows] = idxhit[j];
1488b1683b59SHong Zhang         nrows++;
1489b1683b59SHong Zhang       }
1490b1683b59SHong Zhang     }
1491e88777f2SHong Zhang     den2sp_i += nrows;
1492077f23c2SHong Zhang 
1493b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1494f99a636bSHong Zhang     rows_i += nrows;
1495b1683b59SHong Zhang   }
14960298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1497b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1498b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1499d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1500b28d80bdSHong Zhang 
1501d2853540SHong Zhang   c->colorforrow = colorforrow;
1502d2853540SHong Zhang   c->rows        = rows;
1503f99a636bSHong Zhang   c->den2sp      = den2sp;
1504d2853540SHong Zhang   c->colorforcol = colorforcol;
1505d2853540SHong Zhang   c->columns     = columns;
15062205254eSKarl Rupp 
1507f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1508b1683b59SHong Zhang   PetscFunctionReturn(0);
1509b1683b59SHong Zhang }
1510