xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision cf742fcc470053b8746f003f5618f6b40d3ea504)
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>
10af0996ceSBarry Smith #include <petsc/private/isimpl.h>
1107475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h>
127bab7c10SHong Zhang 
1358cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
14cd093f1dSJed Brown 
155e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
165e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
175e5acdf2Sstefano_zampini #endif
185e5acdf2Sstefano_zampini 
19150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2038baddfdSBarry Smith {
21dfbe8321SBarry Smith   PetscErrorCode ierr;
225e5acdf2Sstefano_zampini #if !defined(PETSC_HAVE_HYPRE)
236540a9cdSHong Zhang   const char     *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
245e5acdf2Sstefano_zampini   PetscInt       nalg = 6;
255e5acdf2Sstefano_zampini #else
265e5acdf2Sstefano_zampini   const char     *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","hypre"};
275e5acdf2Sstefano_zampini   PetscInt       nalg = 7;
285e5acdf2Sstefano_zampini #endif
296540a9cdSHong Zhang   PetscInt       alg = 0; /* set default algorithm */
305c66b693SKris Buschelman 
315c66b693SKris Buschelman   PetscFunctionBegin;
3226be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
33152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
345e5acdf2Sstefano_zampini     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
35d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
363ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
376540a9cdSHong Zhang     switch (alg) {
386540a9cdSHong Zhang     case 1:
39aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
406540a9cdSHong Zhang       break;
416540a9cdSHong Zhang     case 2:
426540a9cdSHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
436540a9cdSHong Zhang       break;
446540a9cdSHong Zhang     case 3:
450ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
466540a9cdSHong Zhang       break;
476540a9cdSHong Zhang     case 4:
488a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
496540a9cdSHong Zhang       break;
506540a9cdSHong Zhang     case 5:
5158cf0668SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
526540a9cdSHong Zhang       break;
535e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
545e5acdf2Sstefano_zampini     case 6:
555e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
565e5acdf2Sstefano_zampini       break;
575e5acdf2Sstefano_zampini #endif
586540a9cdSHong Zhang     default:
5926be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
606540a9cdSHong Zhang      break;
6125023028SHong Zhang     }
623ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
6326be0446SHong Zhang   }
645c913ed7SHong Zhang 
653ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
6601e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
673ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
685c66b693SKris Buschelman   PetscFunctionReturn(0);
695c66b693SKris Buschelman }
701c24bd37SHong Zhang 
7158cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
72b561aa9dSHong Zhang {
73b561aa9dSHong Zhang   PetscErrorCode     ierr;
74b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
75971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
76c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
77b561aa9dSHong Zhang   PetscReal          afill;
78eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
79eca6b86aSHong Zhang   PetscTable         ta;
80fb908031SHong Zhang   PetscBT            lnkbt;
810298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
82b561aa9dSHong Zhang 
83b561aa9dSHong Zhang   PetscFunctionBegin;
84bd958071SHong Zhang   /* Get ci and cj */
85bd958071SHong Zhang   /*---------------*/
86fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
87fb908031SHong Zhang   /* free space for accumulating nonzero column info */
88854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
89fb908031SHong Zhang   ci[0] = 0;
90fb908031SHong Zhang 
91fb908031SHong Zhang   /* create and initialize a linked list */
92c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
93c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
94eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
95eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
96eca6b86aSHong Zhang 
97eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr);
98fb908031SHong Zhang 
99fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
100f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
1012205254eSKarl Rupp 
102fb908031SHong Zhang   current_space = free_space;
103fb908031SHong Zhang 
104fb908031SHong Zhang   /* Determine ci and cj */
105fb908031SHong Zhang   for (i=0; i<am; i++) {
106fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
107fb908031SHong Zhang     aj   = a->j + ai[i];
108fb908031SHong Zhang     for (j=0; j<anzi; j++) {
109fb908031SHong Zhang       brow = aj[j];
110fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
111fb908031SHong Zhang       bj   = b->j + bi[brow];
112fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
113fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
114fb908031SHong Zhang     }
115fb908031SHong Zhang     cnzi = lnk[0];
116fb908031SHong Zhang 
117fb908031SHong Zhang     /* If free space is not available, make more free space */
118fb908031SHong Zhang     /* Double the amount of total space in the list */
119fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
120f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
121fb908031SHong Zhang       ndouble++;
122fb908031SHong Zhang     }
123fb908031SHong Zhang 
124fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
125fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1262205254eSKarl Rupp 
127fb908031SHong Zhang     current_space->array           += cnzi;
128fb908031SHong Zhang     current_space->local_used      += cnzi;
129fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1302205254eSKarl Rupp 
131fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
132fb908031SHong Zhang   }
133fb908031SHong Zhang 
134fb908031SHong Zhang   /* Column indices are in the list of free space */
135fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
136fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
137854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
138fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
139fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
140b561aa9dSHong Zhang 
14126be0446SHong Zhang   /* put together the new symbolic matrix */
142ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
14333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
14458c24d83SHong Zhang 
14558c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
14658c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
14758c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
148aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
149e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
15058c24d83SHong Zhang   c->nonew                  = 0;
151002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1520b7e3e3dSHong Zhang 
1537212da7cSHong Zhang   /* set MatInfo */
1547212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
155f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1567212da7cSHong Zhang   c->maxnz                     = ci[am];
1577212da7cSHong Zhang   c->nz                        = ci[am];
158fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1597212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1607212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1617212da7cSHong Zhang 
1627212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1637212da7cSHong Zhang   if (ci[am]) {
16457622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
16557622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
166f2b054eeSHong Zhang   } else {
167f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
168be0fcf8dSHong Zhang   }
169f2b054eeSHong Zhang #endif
17058c24d83SHong Zhang   PetscFunctionReturn(0);
17158c24d83SHong Zhang }
172d50806bdSBarry Smith 
173dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
174d50806bdSBarry Smith {
175dfbe8321SBarry Smith   PetscErrorCode ierr;
176d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
177d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
178d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
179d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
18038baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
181c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
182519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
183aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
184aa1aec99SHong Zhang   PetscScalar    *ab_dense;
185d50806bdSBarry Smith 
186d50806bdSBarry Smith   PetscFunctionBegin;
187aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
188854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
189aa1aec99SHong Zhang     c->a      = ca;
190aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
191aa1aec99SHong Zhang   } else {
192aa1aec99SHong Zhang     ca        = c->a;
193d90d86d0SMatthew G. Knepley   }
194d90d86d0SMatthew G. Knepley   if (!c->matmult_abdense) {
1951795a4d1SJed Brown     ierr = PetscCalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr);
196d90d86d0SMatthew G. Knepley     c->matmult_abdense = ab_dense;
197d90d86d0SMatthew G. Knepley   } else {
198aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
199aa1aec99SHong Zhang   }
200aa1aec99SHong Zhang 
201c124e916SHong Zhang   /* clean old values in C */
202c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
203d50806bdSBarry Smith   /* Traverse A row-wise. */
204d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
205d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
206d50806bdSBarry Smith   for (i=0; i<am; i++) {
207d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
208d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
209519eb980SHong Zhang       brow = aj[j];
210d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
211d50806bdSBarry Smith       bjj  = bj + bi[brow];
212d50806bdSBarry Smith       baj  = ba + bi[brow];
213519eb980SHong Zhang       /* perform dense axpy */
21436ec6d2dSHong Zhang       valtmp = aa[j];
215519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
21636ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
217519eb980SHong Zhang       }
218519eb980SHong Zhang       flops += 2*bnzi;
219519eb980SHong Zhang     }
220c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
221c58ca2e3SHong Zhang 
222c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
223519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
224519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
225519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
226519eb980SHong Zhang     }
227519eb980SHong Zhang     flops += cnzi;
228519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
229519eb980SHong Zhang   }
230c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
231c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
232c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
233c58ca2e3SHong Zhang   PetscFunctionReturn(0);
234c58ca2e3SHong Zhang }
235c58ca2e3SHong Zhang 
23625023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
237c58ca2e3SHong Zhang {
238c58ca2e3SHong Zhang   PetscErrorCode ierr;
239c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
240c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
241c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
242c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
243c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
244c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
245c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
24636ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
247c58ca2e3SHong Zhang   PetscInt       nextb;
248c58ca2e3SHong Zhang 
249c58ca2e3SHong Zhang   PetscFunctionBegin;
250*cf742fccSHong Zhang   if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
251*cf742fccSHong Zhang     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
252*cf742fccSHong Zhang     c->a      = ca;
253*cf742fccSHong Zhang     c->free_a = PETSC_TRUE;
254*cf742fccSHong Zhang   }
255*cf742fccSHong Zhang 
256c58ca2e3SHong Zhang   /* clean old values in C */
257c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
258c58ca2e3SHong Zhang   /* Traverse A row-wise. */
259c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
260c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
261519eb980SHong Zhang   for (i=0; i<am; i++) {
262519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
263519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
264519eb980SHong Zhang     for (j=0; j<anzi; j++) {
265519eb980SHong Zhang       brow = aj[j];
266519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
267519eb980SHong Zhang       bjj  = bj + bi[brow];
268519eb980SHong Zhang       baj  = ba + bi[brow];
269519eb980SHong Zhang       /* perform sparse axpy */
27036ec6d2dSHong Zhang       valtmp = aa[j];
271c124e916SHong Zhang       nextb  = 0;
272c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
273c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
27436ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
275c124e916SHong Zhang         }
276d50806bdSBarry Smith       }
277d50806bdSBarry Smith       flops += 2*bnzi;
278d50806bdSBarry Smith     }
279519eb980SHong Zhang     aj += anzi; aa += anzi;
280519eb980SHong Zhang     cj += cnzi; ca += cnzi;
281d50806bdSBarry Smith   }
282c58ca2e3SHong Zhang 
283716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
284716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
285d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
286d50806bdSBarry Smith   PetscFunctionReturn(0);
287d50806bdSBarry Smith }
288bc011b1eSHong Zhang 
2893c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
29025296bd5SBarry Smith {
29125296bd5SBarry Smith   PetscErrorCode     ierr;
29225296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
29325296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2943c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
29525296bd5SBarry Smith   MatScalar          *ca;
29625296bd5SBarry Smith   PetscReal          afill;
297eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
298eca6b86aSHong Zhang   PetscTable         ta;
2990298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
30025296bd5SBarry Smith 
30125296bd5SBarry Smith   PetscFunctionBegin;
3023c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
3033c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
3043c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
305854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
3063c50cad2SHong Zhang   ci[0] = 0;
3073c50cad2SHong Zhang 
3083c50cad2SHong Zhang   /* create and initialize a linked list */
309c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
310c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
311eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
312eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
313eca6b86aSHong Zhang 
314eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr);
3153c50cad2SHong Zhang 
3163c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
317f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
3183c50cad2SHong Zhang   current_space = free_space;
3193c50cad2SHong Zhang 
3203c50cad2SHong Zhang   /* Determine ci and cj */
3213c50cad2SHong Zhang   for (i=0; i<am; i++) {
3223c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3233c50cad2SHong Zhang     aj   = a->j + ai[i];
3243c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3253c50cad2SHong Zhang       brow = aj[j];
3263c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3273c50cad2SHong Zhang       bj   = b->j + bi[brow];
3283c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3293c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3303c50cad2SHong Zhang     }
3313c50cad2SHong Zhang     cnzi = lnk[1];
3323c50cad2SHong Zhang 
3333c50cad2SHong Zhang     /* If free space is not available, make more free space */
3343c50cad2SHong Zhang     /* Double the amount of total space in the list */
3353c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
336f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
3373c50cad2SHong Zhang       ndouble++;
3383c50cad2SHong Zhang     }
3393c50cad2SHong Zhang 
3403c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3413c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3422205254eSKarl Rupp 
3433c50cad2SHong Zhang     current_space->array           += cnzi;
3443c50cad2SHong Zhang     current_space->local_used      += cnzi;
3453c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3462205254eSKarl Rupp 
3473c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3483c50cad2SHong Zhang   }
3493c50cad2SHong Zhang 
3503c50cad2SHong Zhang   /* Column indices are in the list of free space */
3513c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3523c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
353854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
3543c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3553c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
35625296bd5SBarry Smith 
35725296bd5SBarry Smith   /* Allocate space for ca */
358854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
35925296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
36025296bd5SBarry Smith 
36125296bd5SBarry Smith   /* put together the new symbolic matrix */
362ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
36333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
36425296bd5SBarry Smith 
36525296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
36625296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
36725296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
36825296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
36925296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
37025296bd5SBarry Smith   c->nonew   = 0;
3712205254eSKarl Rupp 
37225296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
37325296bd5SBarry Smith 
37425296bd5SBarry Smith   /* set MatInfo */
37525296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
37625296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
37725296bd5SBarry Smith   c->maxnz                     = ci[am];
37825296bd5SBarry Smith   c->nz                        = ci[am];
3793c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
38025296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
38125296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
38225296bd5SBarry Smith 
38325296bd5SBarry Smith #if defined(PETSC_USE_INFO)
38425296bd5SBarry Smith   if (ci[am]) {
38557622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
38657622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
38725296bd5SBarry Smith   } else {
38825296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
38925296bd5SBarry Smith   }
39025296bd5SBarry Smith #endif
39125296bd5SBarry Smith   PetscFunctionReturn(0);
39225296bd5SBarry Smith }
39325296bd5SBarry Smith 
39425296bd5SBarry Smith 
39525023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
396e9e4536cSHong Zhang {
397e9e4536cSHong Zhang   PetscErrorCode     ierr;
398e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
399bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
40025c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
401e9e4536cSHong Zhang   MatScalar          *ca;
402e9e4536cSHong Zhang   PetscReal          afill;
403eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
404eca6b86aSHong Zhang   PetscTable         ta;
4050298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
406e9e4536cSHong Zhang 
407e9e4536cSHong Zhang   PetscFunctionBegin;
408bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
409bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
410bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
411854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
412bd958071SHong Zhang   ci[0] = 0;
413bd958071SHong Zhang 
414bd958071SHong Zhang   /* create and initialize a linked list */
415c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
416c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
417eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
418eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
419eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
420bd958071SHong Zhang 
421bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
422f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
423bd958071SHong Zhang   current_space = free_space;
424bd958071SHong Zhang 
425bd958071SHong Zhang   /* Determine ci and cj */
426bd958071SHong Zhang   for (i=0; i<am; i++) {
427bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
428bd958071SHong Zhang     aj   = a->j + ai[i];
429bd958071SHong Zhang     for (j=0; j<anzi; j++) {
430bd958071SHong Zhang       brow = aj[j];
431bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
432bd958071SHong Zhang       bj   = b->j + bi[brow];
433bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
434bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
435bd958071SHong Zhang     }
436bd958071SHong Zhang     cnzi = lnk[0];
437bd958071SHong Zhang 
438bd958071SHong Zhang     /* If free space is not available, make more free space */
439bd958071SHong Zhang     /* Double the amount of total space in the list */
440bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
441f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
442bd958071SHong Zhang       ndouble++;
443bd958071SHong Zhang     }
444bd958071SHong Zhang 
445bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
446bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4472205254eSKarl Rupp 
448bd958071SHong Zhang     current_space->array           += cnzi;
449bd958071SHong Zhang     current_space->local_used      += cnzi;
450bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4512205254eSKarl Rupp 
452bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
453bd958071SHong Zhang   }
454bd958071SHong Zhang 
455bd958071SHong Zhang   /* Column indices are in the list of free space */
456bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
457bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
458854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
459bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
460bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
461e9e4536cSHong Zhang 
462e9e4536cSHong Zhang   /* Allocate space for ca */
463bd958071SHong Zhang   /*-----------------------*/
464854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
465e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
466e9e4536cSHong Zhang 
467e9e4536cSHong Zhang   /* put together the new symbolic matrix */
468ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
46933d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
470e9e4536cSHong Zhang 
471e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
472e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
473e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
474e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
475e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
476e9e4536cSHong Zhang   c->nonew   = 0;
4772205254eSKarl Rupp 
47825023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
479e9e4536cSHong Zhang 
480e9e4536cSHong Zhang   /* set MatInfo */
481e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
482e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
483e9e4536cSHong Zhang   c->maxnz                     = ci[am];
484e9e4536cSHong Zhang   c->nz                        = ci[am];
485bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
486e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
487e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
488e9e4536cSHong Zhang 
489e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
490e9e4536cSHong Zhang   if (ci[am]) {
49157622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
49257622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
493e9e4536cSHong Zhang   } else {
494e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
495e9e4536cSHong Zhang   }
496e9e4536cSHong Zhang #endif
497e9e4536cSHong Zhang   PetscFunctionReturn(0);
498e9e4536cSHong Zhang }
499e9e4536cSHong Zhang 
5000ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
5010ced3a2bSJed Brown {
5020ced3a2bSJed Brown   PetscErrorCode     ierr;
5030ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5040ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
5050ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
5060ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
5070ced3a2bSJed Brown   PetscReal          afill;
5080ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
5090298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
5100ced3a2bSJed Brown   PetscHeap          h;
5110ced3a2bSJed Brown 
5120ced3a2bSJed Brown   PetscFunctionBegin;
513cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
5140ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
5150ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
516854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
5170ced3a2bSJed Brown   ci[0] = 0;
5180ced3a2bSJed Brown 
5190ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
520f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
5210ced3a2bSJed Brown   current_space = free_space;
5220ced3a2bSJed Brown 
5230ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
524785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
5250ced3a2bSJed Brown 
5260ced3a2bSJed Brown   /* Determine ci and cj */
5270ced3a2bSJed Brown   for (i=0; i<am; i++) {
5280ced3a2bSJed 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 */
5290ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5300ced3a2bSJed Brown     ci[i+1] = ci[i];
5310ced3a2bSJed Brown     /* Populate the min heap */
5320ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5330ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5340ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5350ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5360ced3a2bSJed Brown       }
5370ced3a2bSJed Brown     }
5380ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5390ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5400ced3a2bSJed Brown     while (j >= 0) {
5410ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
542f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
5430ced3a2bSJed Brown         ndouble++;
5440ced3a2bSJed Brown       }
5450ced3a2bSJed Brown       *(current_space->array++) = col;
5460ced3a2bSJed Brown       current_space->local_used++;
5470ced3a2bSJed Brown       current_space->local_remaining--;
5480ced3a2bSJed Brown       ci[i+1]++;
5490ced3a2bSJed Brown 
5500ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5510ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5520ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5530ced3a2bSJed Brown         PetscInt j2,col2;
5540ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5550ced3a2bSJed Brown         if (col2 != col) break;
5560ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5570ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5580ced3a2bSJed Brown       }
5590ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5600ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5610ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5620ced3a2bSJed Brown     }
5630ced3a2bSJed Brown   }
5640ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5650ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5660ced3a2bSJed Brown 
5670ced3a2bSJed Brown   /* Column indices are in the list of free space */
5680ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5690ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
570785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
5710ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5720ced3a2bSJed Brown 
5730ced3a2bSJed Brown   /* put together the new symbolic matrix */
574ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
57533d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
5760ced3a2bSJed Brown 
5770ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5780ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5790ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5800ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5810ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5820ced3a2bSJed Brown   c->nonew   = 0;
58326fbe8dcSKarl Rupp 
58489d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5850ced3a2bSJed Brown 
5860ced3a2bSJed Brown   /* set MatInfo */
5870ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5880ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5890ced3a2bSJed Brown   c->maxnz                     = ci[am];
5900ced3a2bSJed Brown   c->nz                        = ci[am];
5910ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5920ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5930ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5940ced3a2bSJed Brown 
5950ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5960ced3a2bSJed Brown   if (ci[am]) {
59757622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
59857622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
5990ced3a2bSJed Brown   } else {
6000ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6010ced3a2bSJed Brown   }
6020ced3a2bSJed Brown #endif
6030ced3a2bSJed Brown   PetscFunctionReturn(0);
6040ced3a2bSJed Brown }
605e9e4536cSHong Zhang 
6068a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
6078a07c6f1SJed Brown {
6088a07c6f1SJed Brown   PetscErrorCode     ierr;
6098a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6108a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
6118a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6128a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6138a07c6f1SJed Brown   PetscReal          afill;
6148a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6150298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6168a07c6f1SJed Brown   PetscHeap          h;
6178a07c6f1SJed Brown   PetscBT            bt;
6188a07c6f1SJed Brown 
6198a07c6f1SJed Brown   PetscFunctionBegin;
620cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6218a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6228a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
623854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
6248a07c6f1SJed Brown   ci[0] = 0;
6258a07c6f1SJed Brown 
6268a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
627f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
6282205254eSKarl Rupp 
6298a07c6f1SJed Brown   current_space = free_space;
6308a07c6f1SJed Brown 
6318a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
632785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
6338a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6348a07c6f1SJed Brown 
6358a07c6f1SJed Brown   /* Determine ci and cj */
6368a07c6f1SJed Brown   for (i=0; i<am; i++) {
6378a07c6f1SJed 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 */
6388a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6398a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6408a07c6f1SJed Brown     ci[i+1] = ci[i];
6418a07c6f1SJed Brown     /* Populate the min heap */
6428a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6438a07c6f1SJed Brown       PetscInt brow = acol[j];
6448a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6458a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6468a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6478a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6488a07c6f1SJed Brown           bb[j]++;
6498a07c6f1SJed Brown           break;
6508a07c6f1SJed Brown         }
6518a07c6f1SJed Brown       }
6528a07c6f1SJed Brown     }
6538a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6548a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6558a07c6f1SJed Brown     while (j >= 0) {
6568a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6570298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
658f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
6598a07c6f1SJed Brown         ndouble++;
6608a07c6f1SJed Brown       }
6618a07c6f1SJed Brown       *(current_space->array++) = col;
6628a07c6f1SJed Brown       current_space->local_used++;
6638a07c6f1SJed Brown       current_space->local_remaining--;
6648a07c6f1SJed Brown       ci[i+1]++;
6658a07c6f1SJed Brown 
6668a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6678a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6688a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6698a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6708a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6718a07c6f1SJed Brown           bb[j]++;
6728a07c6f1SJed Brown           break;
6738a07c6f1SJed Brown         }
6748a07c6f1SJed Brown       }
6758a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6768a07c6f1SJed Brown     }
6778a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6788a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6798a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6808a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6818a07c6f1SJed Brown     }
6828a07c6f1SJed Brown   }
6838a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6848a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6858a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6868a07c6f1SJed Brown 
6878a07c6f1SJed Brown   /* Column indices are in the list of free space */
6888a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6898a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
690785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
6918a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6928a07c6f1SJed Brown 
6938a07c6f1SJed Brown   /* put together the new symbolic matrix */
694ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
69533d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
6968a07c6f1SJed Brown 
6978a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6988a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6998a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
7008a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
7018a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
7028a07c6f1SJed Brown   c->nonew   = 0;
70326fbe8dcSKarl Rupp 
70489d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
7058a07c6f1SJed Brown 
7068a07c6f1SJed Brown   /* set MatInfo */
7078a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
7088a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
7098a07c6f1SJed Brown   c->maxnz                     = ci[am];
7108a07c6f1SJed Brown   c->nz                        = ci[am];
7118a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7128a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7138a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7148a07c6f1SJed Brown 
7158a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7168a07c6f1SJed Brown   if (ci[am]) {
71757622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
71857622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
7198a07c6f1SJed Brown   } else {
7208a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7218a07c6f1SJed Brown   }
7228a07c6f1SJed Brown #endif
7238a07c6f1SJed Brown   PetscFunctionReturn(0);
7248a07c6f1SJed Brown }
7258a07c6f1SJed Brown 
726cd093f1dSJed Brown /* concatenate unique entries and then sort */
72758cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
728cd093f1dSJed Brown {
729cd093f1dSJed Brown   PetscErrorCode     ierr;
730cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
731cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
732cd093f1dSJed Brown   PetscInt           *ci,*cj;
733cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
734cd093f1dSJed Brown   PetscReal          afill;
735cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
736cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
737cd093f1dSJed Brown   char               *seen;
738cd093f1dSJed Brown 
739cd093f1dSJed Brown   PetscFunctionBegin;
740854ce69bSBarry Smith   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
741cd093f1dSJed Brown   ci[0] = 0;
742cd093f1dSJed Brown 
743cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
744cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
745cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
746785e854fSJed Brown   ierr = PetscMalloc1(bn,&seen);CHKERRQ(ierr);
747cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
748cd093f1dSJed Brown 
749cd093f1dSJed Brown   /* Determine ci and cj */
750cd093f1dSJed Brown   for (i=0; i<am; i++) {
751cd093f1dSJed 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 */
752cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
753cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
754cd093f1dSJed Brown     /* Pack segrow */
755cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
756cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
757cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
758cd093f1dSJed Brown         PetscInt bcol = bj[k];
759cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
760cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
761cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
762cd093f1dSJed Brown           *slot = bcol;
763cd093f1dSJed Brown           seen[bcol] = 1;
764cd093f1dSJed Brown           packlen++;
765cd093f1dSJed Brown         }
766cd093f1dSJed Brown       }
767cd093f1dSJed Brown     }
768cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
769cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
770cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
771cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
772cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
773cd093f1dSJed Brown   }
774cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
775cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
776cd093f1dSJed Brown 
777cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
778cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
779cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
780cd093f1dSJed Brown 
781cd093f1dSJed Brown   /* put together the new symbolic matrix */
782cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
78333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
784cd093f1dSJed Brown 
785cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
786cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
787cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
788cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
789cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
790cd093f1dSJed Brown   c->nonew   = 0;
791cd093f1dSJed Brown 
792cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
793cd093f1dSJed Brown 
794cd093f1dSJed Brown   /* set MatInfo */
795cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
796cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
797cd093f1dSJed Brown   c->maxnz                     = ci[am];
798cd093f1dSJed Brown   c->nz                        = ci[am];
799cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
800cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
801cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
802cd093f1dSJed Brown 
803cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
804cd093f1dSJed Brown   if (ci[am]) {
80557622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
80657622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
807cd093f1dSJed Brown   } else {
808cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
809cd093f1dSJed Brown   }
810cd093f1dSJed Brown #endif
811cd093f1dSJed Brown   PetscFunctionReturn(0);
812cd093f1dSJed Brown }
813cd093f1dSJed Brown 
814d2853540SHong Zhang /* This routine is not used. Should be removed! */
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 
8312128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8322128a86cSHong Zhang {
8332128a86cSHong Zhang   PetscErrorCode      ierr;
8344c7df5ccSHong Zhang   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
83540192850SHong Zhang   Mat_MatMatTransMult *abt=a->abt;
8362128a86cSHong Zhang 
8372128a86cSHong Zhang   PetscFunctionBegin;
83840192850SHong Zhang   ierr = (abt->destroy)(A);CHKERRQ(ierr);
83940192850SHong Zhang   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
84040192850SHong Zhang   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
84140192850SHong Zhang   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
84240192850SHong Zhang   ierr = PetscFree(abt);CHKERRQ(ierr);
8432128a86cSHong Zhang   PetscFunctionReturn(0);
8442128a86cSHong Zhang }
8452128a86cSHong Zhang 
8466fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
847bc011b1eSHong Zhang {
8485df89d91SHong Zhang   PetscErrorCode      ierr;
84981d82fe4SHong Zhang   Mat                 Bt;
85081d82fe4SHong Zhang   PetscInt            *bti,*btj;
85140192850SHong Zhang   Mat_MatMatTransMult *abt;
8524c7df5ccSHong Zhang   Mat_SeqAIJ          *c;
853d2853540SHong Zhang 
85481d82fe4SHong Zhang   PetscFunctionBegin;
85581d82fe4SHong Zhang   /* create symbolic Bt */
85681d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8570298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
85833d57670SJed Brown   ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
85981d82fe4SHong Zhang 
86081d82fe4SHong Zhang   /* get symbolic C=A*Bt */
86181d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
86281d82fe4SHong Zhang 
8632128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
864b00a9115SJed Brown   ierr   = PetscNew(&abt);CHKERRQ(ierr);
8654c7df5ccSHong Zhang   c      = (Mat_SeqAIJ*)(*C)->data;
86640192850SHong Zhang   c->abt = abt;
8672128a86cSHong Zhang 
86840192850SHong Zhang   abt->usecoloring = PETSC_FALSE;
86940192850SHong Zhang   abt->destroy     = (*C)->ops->destroy;
8702128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8712128a86cSHong Zhang 
872c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr);
87340192850SHong Zhang   if (abt->usecoloring) {
874b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
875b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
876335efc43SPeter Brune     MatColoring          coloring;
8772128a86cSHong Zhang     ISColoring           iscoloring;
8782128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
8794d478ae7SHong Zhang     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
8804d478ae7SHong Zhang     /* inode causes memory problem, don't know why */
8814d478ae7SHong Zhang     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
882e8354b3eSHong Zhang 
883335efc43SPeter Brune     ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr);
884335efc43SPeter Brune     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
885335efc43SPeter Brune     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
886335efc43SPeter Brune     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
887335efc43SPeter Brune     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
888335efc43SPeter Brune     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
889b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
8902205254eSKarl Rupp 
89140192850SHong Zhang     abt->matcoloring = matcoloring;
8922205254eSKarl Rupp 
8932128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
8942128a86cSHong Zhang 
8952128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
8962128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
8972128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
8982128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
8990298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9002205254eSKarl Rupp 
9012128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
90240192850SHong Zhang     abt->Bt_den   = Bt_dense;
9032128a86cSHong Zhang 
9042128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9052128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9062128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9070298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9082205254eSKarl Rupp 
9092128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
91040192850SHong Zhang     abt->ABt_den  = C_dense;
911f94ccd6cSHong Zhang 
912f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9131ce71dffSSatish Balay     {
914f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
915c40ebe3bSHong 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);
9161ce71dffSSatish Balay     }
917f94ccd6cSHong Zhang #endif
9182128a86cSHong Zhang   }
919e99be685SHong Zhang   /* clean up */
920e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
921e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9225df89d91SHong Zhang   PetscFunctionReturn(0);
9235df89d91SHong Zhang }
9245df89d91SHong Zhang 
9256fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9265df89d91SHong Zhang {
9275df89d91SHong Zhang   PetscErrorCode      ierr;
9285df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
929e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
93081d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9315df89d91SHong Zhang   PetscLogDouble      flops=0.0;
932aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
93340192850SHong Zhang   Mat_MatMatTransMult *abt = c->abt;
9345df89d91SHong Zhang 
9355df89d91SHong Zhang   PetscFunctionBegin;
93658ed3ceaSHong Zhang   /* clear old values in C */
93758ed3ceaSHong Zhang   if (!c->a) {
938854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
93958ed3ceaSHong Zhang     c->a      = ca;
94058ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
94158ed3ceaSHong Zhang   } else {
94258ed3ceaSHong Zhang     ca =  c->a;
94358ed3ceaSHong Zhang   }
94458ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
94558ed3ceaSHong Zhang 
94640192850SHong Zhang   if (abt->usecoloring) {
94740192850SHong Zhang     MatTransposeColoring matcoloring = abt->matcoloring;
94840192850SHong Zhang     Mat                  Bt_dense,C_dense = abt->ABt_den;
949c8db22f6SHong Zhang 
950b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
95140192850SHong Zhang     Bt_dense = abt->Bt_den;
952b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
953c8db22f6SHong Zhang 
954c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9552128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
956c8db22f6SHong Zhang 
9572128a86cSHong Zhang     /* Recover C from C_dense */
958b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
959c8db22f6SHong Zhang     PetscFunctionReturn(0);
960c8db22f6SHong Zhang   }
961c8db22f6SHong Zhang 
9621fa1209cSHong Zhang   for (i=0; i<cm; i++) {
96381d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
96481d82fe4SHong Zhang     acol = aj + ai[i];
9656973769bSHong Zhang     aval = aa + ai[i];
9661fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9671fa1209cSHong Zhang     ccol = cj + ci[i];
9686973769bSHong Zhang     cval = ca + ci[i];
9691fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
97081d82fe4SHong Zhang       brow = ccol[j];
97181d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
97281d82fe4SHong Zhang       bcol = bj + bi[brow];
9736973769bSHong Zhang       bval = ba + bi[brow];
9746973769bSHong Zhang 
97581d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
97681d82fe4SHong Zhang       nexta = 0; nextb = 0;
97781d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
9787b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
97981d82fe4SHong Zhang         if (nexta == anzi) break;
9807b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
98181d82fe4SHong Zhang         if (nextb == bnzj) break;
98281d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
9836973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
98481d82fe4SHong Zhang           nexta++; nextb++;
98581d82fe4SHong Zhang           flops += 2;
9861fa1209cSHong Zhang         }
9871fa1209cSHong Zhang       }
98881d82fe4SHong Zhang     }
98981d82fe4SHong Zhang   }
99081d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99181d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99281d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
9935df89d91SHong Zhang   PetscFunctionReturn(0);
9945df89d91SHong Zhang }
9955df89d91SHong Zhang 
9960adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9970adebc6cSBarry Smith {
9985df89d91SHong Zhang   PetscErrorCode ierr;
9995df89d91SHong Zhang 
10005df89d91SHong Zhang   PetscFunctionBegin;
10015df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
100207706670SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
100375648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
100407706670SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
10055df89d91SHong Zhang   }
100607706670SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
100775648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
100807706670SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
10095df89d91SHong Zhang   PetscFunctionReturn(0);
10105df89d91SHong Zhang }
10115df89d91SHong Zhang 
101275648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10135df89d91SHong Zhang {
1014bc011b1eSHong Zhang   PetscErrorCode ierr;
1015bc011b1eSHong Zhang   Mat            At;
101638baddfdSBarry Smith   PetscInt       *ati,*atj;
1017bc011b1eSHong Zhang 
1018bc011b1eSHong Zhang   PetscFunctionBegin;
1019bc011b1eSHong Zhang   /* create symbolic At */
1020bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10210298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
102233d57670SJed Brown   ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1023bc011b1eSHong Zhang 
1024bc011b1eSHong Zhang   /* get symbolic C=At*B */
1025bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1026bc011b1eSHong Zhang 
1027bc011b1eSHong Zhang   /* clean up */
10286bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1029bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1030bc011b1eSHong Zhang   PetscFunctionReturn(0);
1031bc011b1eSHong Zhang }
1032bc011b1eSHong Zhang 
103375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1034bc011b1eSHong Zhang {
1035bc011b1eSHong Zhang   PetscErrorCode ierr;
10360fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1037d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1038d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1039d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1040aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1041bc011b1eSHong Zhang 
1042bc011b1eSHong Zhang   PetscFunctionBegin;
1043aa1aec99SHong Zhang   if (!c->a) {
1044854ce69bSBarry Smith     ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
10452205254eSKarl Rupp 
1046aa1aec99SHong Zhang     c->a      = ca;
1047aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1048aa1aec99SHong Zhang   } else {
1049aa1aec99SHong Zhang     ca = c->a;
1050aa1aec99SHong Zhang   }
1051bc011b1eSHong Zhang   /* clear old values in C */
1052bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1053bc011b1eSHong Zhang 
1054bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1055bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1056bc011b1eSHong Zhang     bj   = b->j + bi[i];
1057bc011b1eSHong Zhang     ba   = b->a + bi[i];
1058bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1059bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1060bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1061bc011b1eSHong Zhang       nextb = 0;
10620fbc74f4SHong Zhang       crow  = *aj++;
1063bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1064bc011b1eSHong Zhang       caj   = ca + ci[crow];
1065bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1066bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
10670fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
10680fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1069bc011b1eSHong Zhang           nextb++;
1070bc011b1eSHong Zhang         }
1071bc011b1eSHong Zhang       }
1072bc011b1eSHong Zhang       flops += 2*bnzi;
10730fbc74f4SHong Zhang       aa++;
1074bc011b1eSHong Zhang     }
1075bc011b1eSHong Zhang   }
1076bc011b1eSHong Zhang 
1077bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1078bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1079bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1080bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1081bc011b1eSHong Zhang   PetscFunctionReturn(0);
1082bc011b1eSHong Zhang }
10839123193aSHong Zhang 
1084150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10859123193aSHong Zhang {
10869123193aSHong Zhang   PetscErrorCode ierr;
10879123193aSHong Zhang 
10889123193aSHong Zhang   PetscFunctionBegin;
10899123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
10903ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
10919123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
10923ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
10939123193aSHong Zhang   }
10943ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
10959123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
10964614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
10979123193aSHong Zhang   PetscFunctionReturn(0);
10989123193aSHong Zhang }
1099edd81eecSMatthew Knepley 
11009123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11019123193aSHong Zhang {
11029123193aSHong Zhang   PetscErrorCode ierr;
11039123193aSHong Zhang 
11049123193aSHong Zhang   PetscFunctionBegin;
11055a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11062205254eSKarl Rupp 
1107d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11089123193aSHong Zhang   PetscFunctionReturn(0);
11099123193aSHong Zhang }
11109123193aSHong Zhang 
11119123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11129123193aSHong Zhang {
1113f32d5d43SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1114612438f1SStefano Zampini   Mat_SeqDense      *bd = (Mat_SeqDense*)B->data;
11159123193aSHong Zhang   PetscErrorCode    ierr;
1116daf57211SHong Zhang   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1117daf57211SHong Zhang   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1118daf57211SHong Zhang   const PetscInt    *aj;
1119612438f1SStefano Zampini   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1120daf57211SHong Zhang   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
11219123193aSHong Zhang 
11229123193aSHong Zhang   PetscFunctionBegin;
1123f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1124612438f1SStefano 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);
1125e32f2f54SBarry 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);
1126e32f2f54SBarry 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);
1127612438f1SStefano Zampini   b = bd->v;
11288c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1129f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1130730858b9SHong Zhang   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1131f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1132f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1133f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1134f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1135f32d5d43SBarry Smith       aj = a->j + a->i[i];
1136f32d5d43SBarry Smith       aa = a->a + a->i[i];
1137f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1138730858b9SHong Zhang         aatmp = aa[j]; ajtmp = aj[j];
1139730858b9SHong Zhang         r1 += aatmp*b1[ajtmp];
1140730858b9SHong Zhang         r2 += aatmp*b2[ajtmp];
1141730858b9SHong Zhang         r3 += aatmp*b3[ajtmp];
1142730858b9SHong Zhang         r4 += aatmp*b4[ajtmp];
11439123193aSHong Zhang       }
1144730858b9SHong Zhang       c1[i] = r1;
1145730858b9SHong Zhang       c2[i] = r2;
1146730858b9SHong Zhang       c3[i] = r3;
1147730858b9SHong Zhang       c4[i] = r4;
1148f32d5d43SBarry Smith     }
1149730858b9SHong Zhang     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1150730858b9SHong Zhang     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1151f32d5d43SBarry Smith   }
1152f32d5d43SBarry Smith   for (; col<cn; col++) {   /* over extra columns of C */
1153f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1154f32d5d43SBarry Smith       r1 = 0.0;
1155f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1156f32d5d43SBarry Smith       aj = a->j + a->i[i];
1157f32d5d43SBarry Smith       aa = a->a + a->i[i];
1158f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1159730858b9SHong Zhang         r1 += aa[j]*b1[aj[j]];
1160f32d5d43SBarry Smith       }
1161730858b9SHong Zhang       c1[i] = r1;
1162f32d5d43SBarry Smith     }
1163f32d5d43SBarry Smith     b1 += bm;
1164730858b9SHong Zhang     c1 += am;
1165f32d5d43SBarry Smith   }
1166dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
11678c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1168f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1169f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1170f32d5d43SBarry Smith   PetscFunctionReturn(0);
1171f32d5d43SBarry Smith }
1172f32d5d43SBarry Smith 
1173f32d5d43SBarry Smith /*
11744324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1175f32d5d43SBarry Smith */
1176f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1177f32d5d43SBarry Smith {
1178f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
117990f5ea3eSStefano Zampini   Mat_SeqDense   *bd = (Mat_SeqDense*)B->data;
1180f32d5d43SBarry Smith   PetscErrorCode ierr;
1181dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1182dd6ea824SBarry Smith   MatScalar      *aa;
118390f5ea3eSStefano Zampini   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
11844324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1185f32d5d43SBarry Smith 
1186f32d5d43SBarry Smith   PetscFunctionBegin;
1187f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
118890f5ea3eSStefano Zampini   b = bd->v;
11898c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1190f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
11914324174eSBarry Smith 
11924324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
11934324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
11944324174eSBarry Smith       colam = col*am;
11954324174eSBarry Smith       arm   = a->compressedrow.nrows;
11964324174eSBarry Smith       ii    = a->compressedrow.i;
11974324174eSBarry Smith       ridx  = a->compressedrow.rindex;
11984324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
11994324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12004324174eSBarry Smith         n  = ii[i+1] - ii[i];
12014324174eSBarry Smith         aj = a->j + ii[i];
12024324174eSBarry Smith         aa = a->a + ii[i];
12034324174eSBarry Smith         for (j=0; j<n; j++) {
12044324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12054324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12064324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12074324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12084324174eSBarry Smith         }
12094324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12104324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12114324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12124324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12134324174eSBarry Smith       }
12144324174eSBarry Smith       b1 += bm4;
12154324174eSBarry Smith       b2 += bm4;
12164324174eSBarry Smith       b3 += bm4;
12174324174eSBarry Smith       b4 += bm4;
12184324174eSBarry Smith     }
12194324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12204324174eSBarry Smith       colam = col*am;
12214324174eSBarry Smith       arm   = a->compressedrow.nrows;
12224324174eSBarry Smith       ii    = a->compressedrow.i;
12234324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12244324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12254324174eSBarry Smith         r1 = 0.0;
12264324174eSBarry Smith         n  = ii[i+1] - ii[i];
12274324174eSBarry Smith         aj = a->j + ii[i];
12284324174eSBarry Smith         aa = a->a + ii[i];
12294324174eSBarry Smith 
12304324174eSBarry Smith         for (j=0; j<n; j++) {
12314324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
12324324174eSBarry Smith         }
1233a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
12344324174eSBarry Smith       }
12354324174eSBarry Smith       b1 += bm;
12364324174eSBarry Smith     }
12374324174eSBarry Smith   } else {
1238f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1239f32d5d43SBarry Smith       colam = col*am;
1240f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1241f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1242f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1243f32d5d43SBarry Smith         aj = a->j + a->i[i];
1244f32d5d43SBarry Smith         aa = a->a + a->i[i];
1245f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1246f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1247f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1248f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1249f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1250f32d5d43SBarry Smith         }
1251f32d5d43SBarry Smith         c[colam + i]       += r1;
1252f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1253f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1254f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1255f32d5d43SBarry Smith       }
1256f32d5d43SBarry Smith       b1 += bm4;
1257f32d5d43SBarry Smith       b2 += bm4;
1258f32d5d43SBarry Smith       b3 += bm4;
1259f32d5d43SBarry Smith       b4 += bm4;
1260f32d5d43SBarry Smith     }
1261f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1262a2ea699eSBarry Smith       colam = col*am;
1263f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1264f32d5d43SBarry Smith         r1 = 0.0;
1265f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1266f32d5d43SBarry Smith         aj = a->j + a->i[i];
1267f32d5d43SBarry Smith         aa = a->a + a->i[i];
1268f32d5d43SBarry Smith 
1269f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1270f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1271f32d5d43SBarry Smith         }
1272a2ea699eSBarry Smith         c[colam + i] += r1;
1273f32d5d43SBarry Smith       }
1274f32d5d43SBarry Smith       b1 += bm;
1275f32d5d43SBarry Smith     }
12764324174eSBarry Smith   }
1277dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
12788c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
12799123193aSHong Zhang   PetscFunctionReturn(0);
12809123193aSHong Zhang }
1281b1683b59SHong Zhang 
1282b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1283c8db22f6SHong Zhang {
1284c8db22f6SHong Zhang   PetscErrorCode ierr;
12852f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
12862f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
12872f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
12882f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
12892f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
12902f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1291c8db22f6SHong Zhang 
1292c8db22f6SHong Zhang   PetscFunctionBegin;
12932f5f1f90SHong Zhang   btval_den=btdense->v;
12942f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
12952f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
12962f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
12972f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1298d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
12992f5f1f90SHong Zhang       btcol = bj + bi[col];
13002f5f1f90SHong Zhang       btval = ba + bi[col];
13012f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1302d2853540SHong Zhang       for (j=0; j<anz; j++) {
13032f5f1f90SHong Zhang         brow            = btcol[j];
13042f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1305c8db22f6SHong Zhang       }
1306c8db22f6SHong Zhang     }
13072f5f1f90SHong Zhang     btval_den += m;
1308c8db22f6SHong Zhang   }
1309c8db22f6SHong Zhang   PetscFunctionReturn(0);
1310c8db22f6SHong Zhang }
1311c8db22f6SHong Zhang 
1312b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13138972f759SHong Zhang {
13148972f759SHong Zhang   PetscErrorCode ierr;
1315b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
1316077f23c2SHong Zhang   PetscScalar    *ca_den,*ca_den_ptr,*ca=csp->a;
1317f99a636bSHong Zhang   PetscInt       k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1318e88777f2SHong Zhang   PetscInt       brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1319077f23c2SHong Zhang   PetscInt       nrows,*row,*idx;
1320077f23c2SHong Zhang   PetscInt       *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
13218972f759SHong Zhang 
13228972f759SHong Zhang   PetscFunctionBegin;
1323a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1324f99a636bSHong Zhang 
1325077f23c2SHong Zhang   if (brows > 0) {
1326077f23c2SHong Zhang     PetscInt *lstart,row_end,row_start;
1327980ae229SHong Zhang     lstart = matcoloring->lstart;
1328eeb4fd02SHong Zhang     ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr);
1329eeb4fd02SHong Zhang 
1330077f23c2SHong Zhang     row_end = brows;
1331eeb4fd02SHong Zhang     if (row_end > m) row_end = m;
1332077f23c2SHong Zhang     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1333077f23c2SHong Zhang       ca_den_ptr = ca_den;
1334980ae229SHong Zhang       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1335eeb4fd02SHong Zhang         nrows = matcoloring->nrows[k];
1336eeb4fd02SHong Zhang         row   = rows  + colorforrow[k];
1337eeb4fd02SHong Zhang         idx   = den2sp + colorforrow[k];
1338eeb4fd02SHong Zhang         for (l=lstart[k]; l<nrows; l++) {
1339eeb4fd02SHong Zhang           if (row[l] >= row_end) {
1340eeb4fd02SHong Zhang             lstart[k] = l;
1341eeb4fd02SHong Zhang             break;
1342eeb4fd02SHong Zhang           } else {
1343077f23c2SHong Zhang             ca[idx[l]] = ca_den_ptr[row[l]];
1344eeb4fd02SHong Zhang           }
1345eeb4fd02SHong Zhang         }
1346077f23c2SHong Zhang         ca_den_ptr += m;
1347eeb4fd02SHong Zhang       }
1348077f23c2SHong Zhang       row_end += brows;
1349eeb4fd02SHong Zhang       if (row_end > m) row_end = m;
1350eeb4fd02SHong Zhang     }
1351077f23c2SHong Zhang   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1352077f23c2SHong Zhang     ca_den_ptr = ca_den;
1353077f23c2SHong Zhang     for (k=0; k<ncolors; k++) {
1354077f23c2SHong Zhang       nrows = matcoloring->nrows[k];
1355077f23c2SHong Zhang       row   = rows  + colorforrow[k];
1356077f23c2SHong Zhang       idx   = den2sp + colorforrow[k];
1357077f23c2SHong Zhang       for (l=0; l<nrows; l++) {
1358077f23c2SHong Zhang         ca[idx[l]] = ca_den_ptr[row[l]];
1359077f23c2SHong Zhang       }
1360077f23c2SHong Zhang       ca_den_ptr += m;
1361077f23c2SHong Zhang     }
1362f99a636bSHong Zhang   }
1363f99a636bSHong Zhang 
1364a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
1365f99a636bSHong Zhang #if defined(PETSC_USE_INFO)
1366077f23c2SHong Zhang   if (matcoloring->brows > 0) {
1367f99a636bSHong Zhang     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1368e88777f2SHong Zhang   } else {
1369077f23c2SHong Zhang     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1370e88777f2SHong Zhang   }
1371f99a636bSHong Zhang #endif
13728972f759SHong Zhang   PetscFunctionReturn(0);
13738972f759SHong Zhang }
13748972f759SHong Zhang 
1375b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1376b1683b59SHong Zhang {
1377b1683b59SHong Zhang   PetscErrorCode ierr;
1378e88777f2SHong Zhang   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
13791a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1380b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1381b1683b59SHong Zhang   IS             *isa;
1382b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1383e88777f2SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1384e88777f2SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i,brows;
1385e88777f2SHong Zhang   PetscBool      flg;
1386b1683b59SHong Zhang 
1387b1683b59SHong Zhang   PetscFunctionBegin;
1388b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1389e99be685SHong Zhang 
13907ecccc15SHong Zhang   /* bs >1 is not being tested yet! */
1391e88777f2SHong Zhang   Nbs       = mat->cmap->N/bs;
1392b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1393e88777f2SHong Zhang   c->N      = Nbs;
1394e88777f2SHong Zhang   c->m      = c->M;
1395b1683b59SHong Zhang   c->rstart = 0;
1396e88777f2SHong Zhang   c->brows  = 100;
1397b1683b59SHong Zhang 
1398b1683b59SHong Zhang   c->ncolors = nis;
1399dcca6d9dSJed Brown   ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr);
1400854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr);
1401854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr);
1402e88777f2SHong Zhang 
1403e88777f2SHong Zhang   brows = c->brows;
1404c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1405e88777f2SHong Zhang   if (flg) c->brows = brows;
1406eeb4fd02SHong Zhang   if (brows > 0) {
1407854ce69bSBarry Smith     ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr);
1408e88777f2SHong Zhang   }
14092205254eSKarl Rupp 
1410d2853540SHong Zhang   colorforrow[0] = 0;
1411d2853540SHong Zhang   rows_i         = rows;
1412f99a636bSHong Zhang   den2sp_i       = den2sp;
1413b1683b59SHong Zhang 
1414854ce69bSBarry Smith   ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr);
1415854ce69bSBarry Smith   ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr);
14162205254eSKarl Rupp 
1417d2853540SHong Zhang   colorforcol[0] = 0;
1418d2853540SHong Zhang   columns_i      = columns;
1419d2853540SHong Zhang 
1420077f23c2SHong Zhang   /* get column-wise storage of mat */
1421077f23c2SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1422b1683b59SHong Zhang 
1423b28d80bdSHong Zhang   cm   = c->m;
1424854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr);
1425854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr);
1426f99a636bSHong Zhang   for (i=0; i<nis; i++) { /* loop over color */
1427b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1428b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
14292205254eSKarl Rupp 
1430b1683b59SHong Zhang     c->ncolumns[i] = n;
1431b1683b59SHong Zhang     if (n) {
1432d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1433b1683b59SHong Zhang     }
1434d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1435d2853540SHong Zhang     columns_i       += n;
1436b1683b59SHong Zhang 
1437b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1438e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1439e99be685SHong Zhang 
1440b7caf3d6SHong Zhang     for (j=0; j<n; j++) { /* loop over columns*/
1441b1683b59SHong Zhang       col     = is[j];
1442d2853540SHong Zhang       row_idx = cj + ci[col];
1443b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1444b7caf3d6SHong Zhang       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1445e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1446d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1447b1683b59SHong Zhang       }
1448b1683b59SHong Zhang     }
1449b1683b59SHong Zhang     /* count the number of hits */
1450b1683b59SHong Zhang     nrows = 0;
1451e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1452b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1453b1683b59SHong Zhang     }
1454b1683b59SHong Zhang     c->nrows[i]      = nrows;
1455d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1456d2853540SHong Zhang 
1457b1683b59SHong Zhang     nrows = 0;
1458b7caf3d6SHong Zhang     for (j=0; j<cm; j++) { /* loop over rows */
1459b1683b59SHong Zhang       if (rowhit[j]) {
1460d2853540SHong Zhang         rows_i[nrows]   = j;
146112b89a43SHong Zhang         den2sp_i[nrows] = idxhit[j];
1462b1683b59SHong Zhang         nrows++;
1463b1683b59SHong Zhang       }
1464b1683b59SHong Zhang     }
1465e88777f2SHong Zhang     den2sp_i += nrows;
1466077f23c2SHong Zhang 
1467b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1468f99a636bSHong Zhang     rows_i += nrows;
1469b1683b59SHong Zhang   }
14700298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1471b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1472b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1473d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1474b28d80bdSHong Zhang 
1475d2853540SHong Zhang   c->colorforrow = colorforrow;
1476d2853540SHong Zhang   c->rows        = rows;
1477f99a636bSHong Zhang   c->den2sp      = den2sp;
1478d2853540SHong Zhang   c->colorforcol = colorforcol;
1479d2853540SHong Zhang   c->columns     = columns;
14802205254eSKarl Rupp 
1481f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1482b1683b59SHong Zhang   PetscFunctionReturn(0);
1483b1683b59SHong Zhang }
1484