xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 5e5acdf2c7d3246fd94da6d18ffcb8751186fbe1)
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 
16*5e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
17*5e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
18*5e5acdf2Sstefano_zampini #endif
19*5e5acdf2Sstefano_zampini 
206284ec50SHong Zhang #undef __FUNCT__
215c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
22150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2338baddfdSBarry Smith {
24dfbe8321SBarry Smith   PetscErrorCode ierr;
25*5e5acdf2Sstefano_zampini #if !defined(PETSC_HAVE_HYPRE)
266540a9cdSHong Zhang   const char     *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
27*5e5acdf2Sstefano_zampini   PetscInt       nalg = 6;
28*5e5acdf2Sstefano_zampini #else
29*5e5acdf2Sstefano_zampini   const char     *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","hypre"};
30*5e5acdf2Sstefano_zampini   PetscInt       nalg = 7;
31*5e5acdf2Sstefano_zampini #endif
326540a9cdSHong Zhang   PetscInt       alg = 0; /* set default algorithm */
335c66b693SKris Buschelman 
345c66b693SKris Buschelman   PetscFunctionBegin;
3526be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
36152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
37*5e5acdf2Sstefano_zampini     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
38d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
393ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
406540a9cdSHong Zhang     switch (alg) {
416540a9cdSHong Zhang     case 1:
42aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
436540a9cdSHong Zhang       break;
446540a9cdSHong Zhang     case 2:
456540a9cdSHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
466540a9cdSHong Zhang       break;
476540a9cdSHong Zhang     case 3:
480ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
496540a9cdSHong Zhang       break;
506540a9cdSHong Zhang     case 4:
518a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
526540a9cdSHong Zhang       break;
536540a9cdSHong Zhang     case 5:
5458cf0668SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
556540a9cdSHong Zhang       break;
56*5e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
57*5e5acdf2Sstefano_zampini     case 6:
58*5e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
59*5e5acdf2Sstefano_zampini       break;
60*5e5acdf2Sstefano_zampini #endif
616540a9cdSHong Zhang     default:
6226be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
636540a9cdSHong Zhang      break;
6425023028SHong Zhang     }
653ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
6626be0446SHong Zhang   }
675c913ed7SHong Zhang 
683ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
6901e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
703ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
715c66b693SKris Buschelman   PetscFunctionReturn(0);
725c66b693SKris Buschelman }
731c24bd37SHong Zhang 
74e9e4536cSHong Zhang #undef __FUNCT__
7558cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed"
7658cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
77b561aa9dSHong Zhang {
78b561aa9dSHong Zhang   PetscErrorCode     ierr;
79b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
80971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
81c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
82b561aa9dSHong Zhang   PetscReal          afill;
83eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
84eca6b86aSHong Zhang   PetscTable         ta;
85fb908031SHong Zhang   PetscBT            lnkbt;
860298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
87b561aa9dSHong Zhang 
88b561aa9dSHong Zhang   PetscFunctionBegin;
89bd958071SHong Zhang   /* Get ci and cj */
90bd958071SHong Zhang   /*---------------*/
91fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
92fb908031SHong Zhang   /* free space for accumulating nonzero column info */
93854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
94fb908031SHong Zhang   ci[0] = 0;
95fb908031SHong Zhang 
96fb908031SHong Zhang   /* create and initialize a linked list */
97c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
98c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
99eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
100eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
101eca6b86aSHong Zhang 
102eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr);
103fb908031SHong Zhang 
104fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
105f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
1062205254eSKarl Rupp 
107fb908031SHong Zhang   current_space = free_space;
108fb908031SHong Zhang 
109fb908031SHong Zhang   /* Determine ci and cj */
110fb908031SHong Zhang   for (i=0; i<am; i++) {
111fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
112fb908031SHong Zhang     aj   = a->j + ai[i];
113fb908031SHong Zhang     for (j=0; j<anzi; j++) {
114fb908031SHong Zhang       brow = aj[j];
115fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
116fb908031SHong Zhang       bj   = b->j + bi[brow];
117fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
118fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
119fb908031SHong Zhang     }
120fb908031SHong Zhang     cnzi = lnk[0];
121fb908031SHong Zhang 
122fb908031SHong Zhang     /* If free space is not available, make more free space */
123fb908031SHong Zhang     /* Double the amount of total space in the list */
124fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
125f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
126fb908031SHong Zhang       ndouble++;
127fb908031SHong Zhang     }
128fb908031SHong Zhang 
129fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
130fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1312205254eSKarl Rupp 
132fb908031SHong Zhang     current_space->array           += cnzi;
133fb908031SHong Zhang     current_space->local_used      += cnzi;
134fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1352205254eSKarl Rupp 
136fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
137fb908031SHong Zhang   }
138fb908031SHong Zhang 
139fb908031SHong Zhang   /* Column indices are in the list of free space */
140fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
141fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
142854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
143fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
144fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
145b561aa9dSHong Zhang 
14626be0446SHong Zhang   /* put together the new symbolic matrix */
147ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
14833d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
14958c24d83SHong Zhang 
15058c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
15158c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
15258c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
153aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
154e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
15558c24d83SHong Zhang   c->nonew                  = 0;
156002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1570b7e3e3dSHong Zhang 
1587212da7cSHong Zhang   /* set MatInfo */
1597212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
160f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1617212da7cSHong Zhang   c->maxnz                     = ci[am];
1627212da7cSHong Zhang   c->nz                        = ci[am];
163fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1647212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1657212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1667212da7cSHong Zhang 
1677212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1687212da7cSHong Zhang   if (ci[am]) {
16957622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
17057622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
171f2b054eeSHong Zhang   } else {
172f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
173be0fcf8dSHong Zhang   }
174f2b054eeSHong Zhang #endif
17558c24d83SHong Zhang   PetscFunctionReturn(0);
17658c24d83SHong Zhang }
177d50806bdSBarry Smith 
17826be0446SHong Zhang #undef __FUNCT__
17926be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
180dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
181d50806bdSBarry Smith {
182dfbe8321SBarry Smith   PetscErrorCode ierr;
183d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
184d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
185d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
186d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
18738baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
188c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
189519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
190aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
191aa1aec99SHong Zhang   PetscScalar    *ab_dense;
192d50806bdSBarry Smith 
193d50806bdSBarry Smith   PetscFunctionBegin;
194aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
195854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
196aa1aec99SHong Zhang     c->a      = ca;
197aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
198aa1aec99SHong Zhang   } else {
199aa1aec99SHong Zhang     ca        = c->a;
200d90d86d0SMatthew G. Knepley   }
201d90d86d0SMatthew G. Knepley   if (!c->matmult_abdense) {
2021795a4d1SJed Brown     ierr = PetscCalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr);
203d90d86d0SMatthew G. Knepley     c->matmult_abdense = ab_dense;
204d90d86d0SMatthew G. Knepley   } else {
205aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
206aa1aec99SHong Zhang   }
207aa1aec99SHong Zhang 
208c124e916SHong Zhang   /* clean old values in C */
209c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
210d50806bdSBarry Smith   /* Traverse A row-wise. */
211d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
212d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
213d50806bdSBarry Smith   for (i=0; i<am; i++) {
214d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
215d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
216519eb980SHong Zhang       brow = aj[j];
217d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
218d50806bdSBarry Smith       bjj  = bj + bi[brow];
219d50806bdSBarry Smith       baj  = ba + bi[brow];
220519eb980SHong Zhang       /* perform dense axpy */
22136ec6d2dSHong Zhang       valtmp = aa[j];
222519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
22336ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
224519eb980SHong Zhang       }
225519eb980SHong Zhang       flops += 2*bnzi;
226519eb980SHong Zhang     }
227c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
228c58ca2e3SHong Zhang 
229c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
230519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
231519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
232519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
233519eb980SHong Zhang     }
234519eb980SHong Zhang     flops += cnzi;
235519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
236519eb980SHong Zhang   }
237c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
238c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
239c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
240c58ca2e3SHong Zhang   PetscFunctionReturn(0);
241c58ca2e3SHong Zhang }
242c58ca2e3SHong Zhang 
243c58ca2e3SHong Zhang #undef __FUNCT__
24425023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
24525023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
246c58ca2e3SHong Zhang {
247c58ca2e3SHong Zhang   PetscErrorCode ierr;
248c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
249c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
250c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
251c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
252c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
253c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
254c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
25536ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
256c58ca2e3SHong Zhang   PetscInt       nextb;
257c58ca2e3SHong Zhang 
258c58ca2e3SHong Zhang   PetscFunctionBegin;
259c58ca2e3SHong Zhang   /* clean old values in C */
260c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
261c58ca2e3SHong Zhang   /* Traverse A row-wise. */
262c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
263c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
264519eb980SHong Zhang   for (i=0; i<am; i++) {
265519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
266519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
267519eb980SHong Zhang     for (j=0; j<anzi; j++) {
268519eb980SHong Zhang       brow = aj[j];
269519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
270519eb980SHong Zhang       bjj  = bj + bi[brow];
271519eb980SHong Zhang       baj  = ba + bi[brow];
272519eb980SHong Zhang       /* perform sparse axpy */
27336ec6d2dSHong Zhang       valtmp = aa[j];
274c124e916SHong Zhang       nextb  = 0;
275c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
276c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
27736ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
278c124e916SHong Zhang         }
279d50806bdSBarry Smith       }
280d50806bdSBarry Smith       flops += 2*bnzi;
281d50806bdSBarry Smith     }
282519eb980SHong Zhang     aj += anzi; aa += anzi;
283519eb980SHong Zhang     cj += cnzi; ca += cnzi;
284d50806bdSBarry Smith   }
285c58ca2e3SHong Zhang 
286716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
287716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
288d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
289d50806bdSBarry Smith   PetscFunctionReturn(0);
290d50806bdSBarry Smith }
291bc011b1eSHong Zhang 
292e9e4536cSHong Zhang #undef __FUNCT__
2933c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2943c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
29525296bd5SBarry Smith {
29625296bd5SBarry Smith   PetscErrorCode     ierr;
29725296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
29825296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2993c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
30025296bd5SBarry Smith   MatScalar          *ca;
30125296bd5SBarry Smith   PetscReal          afill;
302eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
303eca6b86aSHong Zhang   PetscTable         ta;
3040298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
30525296bd5SBarry Smith 
30625296bd5SBarry Smith   PetscFunctionBegin;
3073c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
3083c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
3093c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
310854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
3113c50cad2SHong Zhang   ci[0] = 0;
3123c50cad2SHong Zhang 
3133c50cad2SHong Zhang   /* create and initialize a linked list */
314c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
315c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
316eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
317eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
318eca6b86aSHong Zhang 
319eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr);
3203c50cad2SHong Zhang 
3213c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
322f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
3233c50cad2SHong Zhang   current_space = free_space;
3243c50cad2SHong Zhang 
3253c50cad2SHong Zhang   /* Determine ci and cj */
3263c50cad2SHong Zhang   for (i=0; i<am; i++) {
3273c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3283c50cad2SHong Zhang     aj   = a->j + ai[i];
3293c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3303c50cad2SHong Zhang       brow = aj[j];
3313c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3323c50cad2SHong Zhang       bj   = b->j + bi[brow];
3333c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3343c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3353c50cad2SHong Zhang     }
3363c50cad2SHong Zhang     cnzi = lnk[1];
3373c50cad2SHong Zhang 
3383c50cad2SHong Zhang     /* If free space is not available, make more free space */
3393c50cad2SHong Zhang     /* Double the amount of total space in the list */
3403c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
341f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
3423c50cad2SHong Zhang       ndouble++;
3433c50cad2SHong Zhang     }
3443c50cad2SHong Zhang 
3453c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3463c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3472205254eSKarl Rupp 
3483c50cad2SHong Zhang     current_space->array           += cnzi;
3493c50cad2SHong Zhang     current_space->local_used      += cnzi;
3503c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3512205254eSKarl Rupp 
3523c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3533c50cad2SHong Zhang   }
3543c50cad2SHong Zhang 
3553c50cad2SHong Zhang   /* Column indices are in the list of free space */
3563c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3573c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
358854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
3593c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3603c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
36125296bd5SBarry Smith 
36225296bd5SBarry Smith   /* Allocate space for ca */
363854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
36425296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
36525296bd5SBarry Smith 
36625296bd5SBarry Smith   /* put together the new symbolic matrix */
367ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
36833d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
36925296bd5SBarry Smith 
37025296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
37125296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
37225296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
37325296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
37425296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
37525296bd5SBarry Smith   c->nonew   = 0;
3762205254eSKarl Rupp 
37725296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
37825296bd5SBarry Smith 
37925296bd5SBarry Smith   /* set MatInfo */
38025296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
38125296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
38225296bd5SBarry Smith   c->maxnz                     = ci[am];
38325296bd5SBarry Smith   c->nz                        = ci[am];
3843c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
38525296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
38625296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
38725296bd5SBarry Smith 
38825296bd5SBarry Smith #if defined(PETSC_USE_INFO)
38925296bd5SBarry Smith   if (ci[am]) {
39057622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
39157622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
39225296bd5SBarry Smith   } else {
39325296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
39425296bd5SBarry Smith   }
39525296bd5SBarry Smith #endif
39625296bd5SBarry Smith   PetscFunctionReturn(0);
39725296bd5SBarry Smith }
39825296bd5SBarry Smith 
39925296bd5SBarry Smith 
40025296bd5SBarry Smith #undef __FUNCT__
40125023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
40225023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
403e9e4536cSHong Zhang {
404e9e4536cSHong Zhang   PetscErrorCode     ierr;
405e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
406bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
40725c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
408e9e4536cSHong Zhang   MatScalar          *ca;
409e9e4536cSHong Zhang   PetscReal          afill;
410eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
411eca6b86aSHong Zhang   PetscTable         ta;
4120298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
413e9e4536cSHong Zhang 
414e9e4536cSHong Zhang   PetscFunctionBegin;
415bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
416bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
417bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
418854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
419bd958071SHong Zhang   ci[0] = 0;
420bd958071SHong Zhang 
421bd958071SHong Zhang   /* create and initialize a linked list */
422c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
423c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
424eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
425eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
426eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
427bd958071SHong Zhang 
428bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
429f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
430bd958071SHong Zhang   current_space = free_space;
431bd958071SHong Zhang 
432bd958071SHong Zhang   /* Determine ci and cj */
433bd958071SHong Zhang   for (i=0; i<am; i++) {
434bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
435bd958071SHong Zhang     aj   = a->j + ai[i];
436bd958071SHong Zhang     for (j=0; j<anzi; j++) {
437bd958071SHong Zhang       brow = aj[j];
438bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
439bd958071SHong Zhang       bj   = b->j + bi[brow];
440bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
441bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
442bd958071SHong Zhang     }
443bd958071SHong Zhang     cnzi = lnk[0];
444bd958071SHong Zhang 
445bd958071SHong Zhang     /* If free space is not available, make more free space */
446bd958071SHong Zhang     /* Double the amount of total space in the list */
447bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
448f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
449bd958071SHong Zhang       ndouble++;
450bd958071SHong Zhang     }
451bd958071SHong Zhang 
452bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
453bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4542205254eSKarl Rupp 
455bd958071SHong Zhang     current_space->array           += cnzi;
456bd958071SHong Zhang     current_space->local_used      += cnzi;
457bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4582205254eSKarl Rupp 
459bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
460bd958071SHong Zhang   }
461bd958071SHong Zhang 
462bd958071SHong Zhang   /* Column indices are in the list of free space */
463bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
464bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
465854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
466bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
467bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
468e9e4536cSHong Zhang 
469e9e4536cSHong Zhang   /* Allocate space for ca */
470bd958071SHong Zhang   /*-----------------------*/
471854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
472e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
473e9e4536cSHong Zhang 
474e9e4536cSHong Zhang   /* put together the new symbolic matrix */
475ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
47633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
477e9e4536cSHong Zhang 
478e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
479e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
480e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
481e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
482e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
483e9e4536cSHong Zhang   c->nonew   = 0;
4842205254eSKarl Rupp 
48525023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
486e9e4536cSHong Zhang 
487e9e4536cSHong Zhang   /* set MatInfo */
488e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
489e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
490e9e4536cSHong Zhang   c->maxnz                     = ci[am];
491e9e4536cSHong Zhang   c->nz                        = ci[am];
492bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
493e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
494e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
495e9e4536cSHong Zhang 
496e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
497e9e4536cSHong Zhang   if (ci[am]) {
49857622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
49957622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
500e9e4536cSHong Zhang   } else {
501e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
502e9e4536cSHong Zhang   }
503e9e4536cSHong Zhang #endif
504e9e4536cSHong Zhang   PetscFunctionReturn(0);
505e9e4536cSHong Zhang }
506e9e4536cSHong Zhang 
5070ced3a2bSJed Brown #undef __FUNCT__
5080ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
5090ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
5100ced3a2bSJed Brown {
5110ced3a2bSJed Brown   PetscErrorCode     ierr;
5120ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5130ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
5140ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
5150ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
5160ced3a2bSJed Brown   PetscReal          afill;
5170ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
5180298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
5190ced3a2bSJed Brown   PetscHeap          h;
5200ced3a2bSJed Brown 
5210ced3a2bSJed Brown   PetscFunctionBegin;
522cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
5230ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
5240ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
525854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
5260ced3a2bSJed Brown   ci[0] = 0;
5270ced3a2bSJed Brown 
5280ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
529f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
5300ced3a2bSJed Brown   current_space = free_space;
5310ced3a2bSJed Brown 
5320ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
533785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
5340ced3a2bSJed Brown 
5350ced3a2bSJed Brown   /* Determine ci and cj */
5360ced3a2bSJed Brown   for (i=0; i<am; i++) {
5370ced3a2bSJed 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 */
5380ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5390ced3a2bSJed Brown     ci[i+1] = ci[i];
5400ced3a2bSJed Brown     /* Populate the min heap */
5410ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5420ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5430ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5440ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5450ced3a2bSJed Brown       }
5460ced3a2bSJed Brown     }
5470ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5480ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5490ced3a2bSJed Brown     while (j >= 0) {
5500ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
551f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
5520ced3a2bSJed Brown         ndouble++;
5530ced3a2bSJed Brown       }
5540ced3a2bSJed Brown       *(current_space->array++) = col;
5550ced3a2bSJed Brown       current_space->local_used++;
5560ced3a2bSJed Brown       current_space->local_remaining--;
5570ced3a2bSJed Brown       ci[i+1]++;
5580ced3a2bSJed Brown 
5590ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5600ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5610ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5620ced3a2bSJed Brown         PetscInt j2,col2;
5630ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5640ced3a2bSJed Brown         if (col2 != col) break;
5650ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5660ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5670ced3a2bSJed Brown       }
5680ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5690ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5700ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5710ced3a2bSJed Brown     }
5720ced3a2bSJed Brown   }
5730ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5740ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5750ced3a2bSJed Brown 
5760ced3a2bSJed Brown   /* Column indices are in the list of free space */
5770ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5780ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
579785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
5800ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5810ced3a2bSJed Brown 
5820ced3a2bSJed Brown   /* put together the new symbolic matrix */
583ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
58433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
5850ced3a2bSJed Brown 
5860ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5870ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5880ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5890ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5900ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5910ced3a2bSJed Brown   c->nonew   = 0;
59226fbe8dcSKarl Rupp 
59389d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5940ced3a2bSJed Brown 
5950ced3a2bSJed Brown   /* set MatInfo */
5960ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5970ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5980ced3a2bSJed Brown   c->maxnz                     = ci[am];
5990ced3a2bSJed Brown   c->nz                        = ci[am];
6000ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
6010ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
6020ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
6030ced3a2bSJed Brown 
6040ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
6050ced3a2bSJed Brown   if (ci[am]) {
60657622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
60757622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
6080ced3a2bSJed Brown   } else {
6090ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6100ced3a2bSJed Brown   }
6110ced3a2bSJed Brown #endif
6120ced3a2bSJed Brown   PetscFunctionReturn(0);
6130ced3a2bSJed Brown }
614e9e4536cSHong Zhang 
6158a07c6f1SJed Brown #undef __FUNCT__
6168a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
6178a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
6188a07c6f1SJed Brown {
6198a07c6f1SJed Brown   PetscErrorCode     ierr;
6208a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6218a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
6228a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6238a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6248a07c6f1SJed Brown   PetscReal          afill;
6258a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6260298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6278a07c6f1SJed Brown   PetscHeap          h;
6288a07c6f1SJed Brown   PetscBT            bt;
6298a07c6f1SJed Brown 
6308a07c6f1SJed Brown   PetscFunctionBegin;
631cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6328a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6338a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
634854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
6358a07c6f1SJed Brown   ci[0] = 0;
6368a07c6f1SJed Brown 
6378a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
638f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
6392205254eSKarl Rupp 
6408a07c6f1SJed Brown   current_space = free_space;
6418a07c6f1SJed Brown 
6428a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
643785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
6448a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6458a07c6f1SJed Brown 
6468a07c6f1SJed Brown   /* Determine ci and cj */
6478a07c6f1SJed Brown   for (i=0; i<am; i++) {
6488a07c6f1SJed 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 */
6498a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6508a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6518a07c6f1SJed Brown     ci[i+1] = ci[i];
6528a07c6f1SJed Brown     /* Populate the min heap */
6538a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6548a07c6f1SJed Brown       PetscInt brow = acol[j];
6558a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6568a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6578a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6588a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6598a07c6f1SJed Brown           bb[j]++;
6608a07c6f1SJed Brown           break;
6618a07c6f1SJed Brown         }
6628a07c6f1SJed Brown       }
6638a07c6f1SJed Brown     }
6648a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6658a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6668a07c6f1SJed Brown     while (j >= 0) {
6678a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6680298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
669f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
6708a07c6f1SJed Brown         ndouble++;
6718a07c6f1SJed Brown       }
6728a07c6f1SJed Brown       *(current_space->array++) = col;
6738a07c6f1SJed Brown       current_space->local_used++;
6748a07c6f1SJed Brown       current_space->local_remaining--;
6758a07c6f1SJed Brown       ci[i+1]++;
6768a07c6f1SJed Brown 
6778a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6788a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6798a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6808a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6818a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6828a07c6f1SJed Brown           bb[j]++;
6838a07c6f1SJed Brown           break;
6848a07c6f1SJed Brown         }
6858a07c6f1SJed Brown       }
6868a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6878a07c6f1SJed Brown     }
6888a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6898a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6908a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6918a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6928a07c6f1SJed Brown     }
6938a07c6f1SJed Brown   }
6948a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6958a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6968a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6978a07c6f1SJed Brown 
6988a07c6f1SJed Brown   /* Column indices are in the list of free space */
6998a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
7008a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
701785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
7028a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
7038a07c6f1SJed Brown 
7048a07c6f1SJed Brown   /* put together the new symbolic matrix */
705ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
70633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
7078a07c6f1SJed Brown 
7088a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
7098a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
7108a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
7118a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
7128a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
7138a07c6f1SJed Brown   c->nonew   = 0;
71426fbe8dcSKarl Rupp 
71589d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
7168a07c6f1SJed Brown 
7178a07c6f1SJed Brown   /* set MatInfo */
7188a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
7198a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
7208a07c6f1SJed Brown   c->maxnz                     = ci[am];
7218a07c6f1SJed Brown   c->nz                        = ci[am];
7228a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7238a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7248a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7258a07c6f1SJed Brown 
7268a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7278a07c6f1SJed Brown   if (ci[am]) {
72857622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
72957622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
7308a07c6f1SJed Brown   } else {
7318a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7328a07c6f1SJed Brown   }
7338a07c6f1SJed Brown #endif
7348a07c6f1SJed Brown   PetscFunctionReturn(0);
7358a07c6f1SJed Brown }
7368a07c6f1SJed Brown 
737cd093f1dSJed Brown #undef __FUNCT__
73858cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
739cd093f1dSJed Brown /* concatenate unique entries and then sort */
74058cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
741cd093f1dSJed Brown {
742cd093f1dSJed Brown   PetscErrorCode     ierr;
743cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
744cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
745cd093f1dSJed Brown   PetscInt           *ci,*cj;
746cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
747cd093f1dSJed Brown   PetscReal          afill;
748cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
749cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
750cd093f1dSJed Brown   char               *seen;
751cd093f1dSJed Brown 
752cd093f1dSJed Brown   PetscFunctionBegin;
753854ce69bSBarry Smith   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
754cd093f1dSJed Brown   ci[0] = 0;
755cd093f1dSJed Brown 
756cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
757cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
758cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
759785e854fSJed Brown   ierr = PetscMalloc1(bn,&seen);CHKERRQ(ierr);
760cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
761cd093f1dSJed Brown 
762cd093f1dSJed Brown   /* Determine ci and cj */
763cd093f1dSJed Brown   for (i=0; i<am; i++) {
764cd093f1dSJed 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 */
765cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
766cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
767cd093f1dSJed Brown     /* Pack segrow */
768cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
769cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
770cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
771cd093f1dSJed Brown         PetscInt bcol = bj[k];
772cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
773cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
774cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
775cd093f1dSJed Brown           *slot = bcol;
776cd093f1dSJed Brown           seen[bcol] = 1;
777cd093f1dSJed Brown           packlen++;
778cd093f1dSJed Brown         }
779cd093f1dSJed Brown       }
780cd093f1dSJed Brown     }
781cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
782cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
783cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
784cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
785cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
786cd093f1dSJed Brown   }
787cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
788cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
789cd093f1dSJed Brown 
790cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
791cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
792cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
793cd093f1dSJed Brown 
794cd093f1dSJed Brown   /* put together the new symbolic matrix */
795cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
79633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
797cd093f1dSJed Brown 
798cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
799cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
800cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
801cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
802cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
803cd093f1dSJed Brown   c->nonew   = 0;
804cd093f1dSJed Brown 
805cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
806cd093f1dSJed Brown 
807cd093f1dSJed Brown   /* set MatInfo */
808cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
809cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
810cd093f1dSJed Brown   c->maxnz                     = ci[am];
811cd093f1dSJed Brown   c->nz                        = ci[am];
812cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
813cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
814cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
815cd093f1dSJed Brown 
816cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
817cd093f1dSJed Brown   if (ci[am]) {
81857622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
81957622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
820cd093f1dSJed Brown   } else {
821cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
822cd093f1dSJed Brown   }
823cd093f1dSJed Brown #endif
824cd093f1dSJed Brown   PetscFunctionReturn(0);
825cd093f1dSJed Brown }
826cd093f1dSJed Brown 
827d2853540SHong Zhang /* This routine is not used. Should be removed! */
828bc011b1eSHong Zhang #undef __FUNCT__
8296fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
8306fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8315df89d91SHong Zhang {
832bc011b1eSHong Zhang   PetscErrorCode ierr;
833bc011b1eSHong Zhang 
834bc011b1eSHong Zhang   PetscFunctionBegin;
835bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
8363ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
8376fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8383ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
839bc011b1eSHong Zhang   }
8403ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
8416fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8423ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
843bc011b1eSHong Zhang   PetscFunctionReturn(0);
844bc011b1eSHong Zhang }
845bc011b1eSHong Zhang 
846bc011b1eSHong Zhang #undef __FUNCT__
8472128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
8482128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8492128a86cSHong Zhang {
8502128a86cSHong Zhang   PetscErrorCode      ierr;
8514c7df5ccSHong Zhang   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
85240192850SHong Zhang   Mat_MatMatTransMult *abt=a->abt;
8532128a86cSHong Zhang 
8542128a86cSHong Zhang   PetscFunctionBegin;
85540192850SHong Zhang   ierr = (abt->destroy)(A);CHKERRQ(ierr);
85640192850SHong Zhang   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
85740192850SHong Zhang   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
85840192850SHong Zhang   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
85940192850SHong Zhang   ierr = PetscFree(abt);CHKERRQ(ierr);
8602128a86cSHong Zhang   PetscFunctionReturn(0);
8612128a86cSHong Zhang }
8622128a86cSHong Zhang 
8632128a86cSHong Zhang #undef __FUNCT__
8646fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
8656fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
866bc011b1eSHong Zhang {
8675df89d91SHong Zhang   PetscErrorCode      ierr;
86881d82fe4SHong Zhang   Mat                 Bt;
86981d82fe4SHong Zhang   PetscInt            *bti,*btj;
87040192850SHong Zhang   Mat_MatMatTransMult *abt;
8714c7df5ccSHong Zhang   Mat_SeqAIJ          *c;
872d2853540SHong Zhang 
87381d82fe4SHong Zhang   PetscFunctionBegin;
87481d82fe4SHong Zhang   /* create symbolic Bt */
87581d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8760298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
87733d57670SJed Brown   ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
87881d82fe4SHong Zhang 
87981d82fe4SHong Zhang   /* get symbolic C=A*Bt */
88081d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
88181d82fe4SHong Zhang 
8822128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
883b00a9115SJed Brown   ierr   = PetscNew(&abt);CHKERRQ(ierr);
8844c7df5ccSHong Zhang   c      = (Mat_SeqAIJ*)(*C)->data;
88540192850SHong Zhang   c->abt = abt;
8862128a86cSHong Zhang 
88740192850SHong Zhang   abt->usecoloring = PETSC_FALSE;
88840192850SHong Zhang   abt->destroy     = (*C)->ops->destroy;
8892128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8902128a86cSHong Zhang 
891c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr);
89240192850SHong Zhang   if (abt->usecoloring) {
893b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
894b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
895335efc43SPeter Brune     MatColoring          coloring;
8962128a86cSHong Zhang     ISColoring           iscoloring;
8972128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
8984d478ae7SHong Zhang     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
8994d478ae7SHong Zhang     /* inode causes memory problem, don't know why */
9004d478ae7SHong Zhang     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
901e8354b3eSHong Zhang 
902335efc43SPeter Brune     ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr);
903335efc43SPeter Brune     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
904335efc43SPeter Brune     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
905335efc43SPeter Brune     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
906335efc43SPeter Brune     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
907335efc43SPeter Brune     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
908b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
9092205254eSKarl Rupp 
91040192850SHong Zhang     abt->matcoloring = matcoloring;
9112205254eSKarl Rupp 
9122128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
9132128a86cSHong Zhang 
9142128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
9152128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
9162128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9172128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
9180298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9192205254eSKarl Rupp 
9202128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
92140192850SHong Zhang     abt->Bt_den   = Bt_dense;
9222128a86cSHong Zhang 
9232128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9242128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9252128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9260298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9272205254eSKarl Rupp 
9282128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
92940192850SHong Zhang     abt->ABt_den  = C_dense;
930f94ccd6cSHong Zhang 
931f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9321ce71dffSSatish Balay     {
933f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
934c40ebe3bSHong 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);
9351ce71dffSSatish Balay     }
936f94ccd6cSHong Zhang #endif
9372128a86cSHong Zhang   }
938e99be685SHong Zhang   /* clean up */
939e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
940e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9415df89d91SHong Zhang   PetscFunctionReturn(0);
9425df89d91SHong Zhang }
9435df89d91SHong Zhang 
9445df89d91SHong Zhang #undef __FUNCT__
9456fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9466fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9475df89d91SHong Zhang {
9485df89d91SHong Zhang   PetscErrorCode      ierr;
9495df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
950e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
95181d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9525df89d91SHong Zhang   PetscLogDouble      flops=0.0;
953aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
95440192850SHong Zhang   Mat_MatMatTransMult *abt = c->abt;
9555df89d91SHong Zhang 
9565df89d91SHong Zhang   PetscFunctionBegin;
95758ed3ceaSHong Zhang   /* clear old values in C */
95858ed3ceaSHong Zhang   if (!c->a) {
959854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
96058ed3ceaSHong Zhang     c->a      = ca;
96158ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
96258ed3ceaSHong Zhang   } else {
96358ed3ceaSHong Zhang     ca =  c->a;
96458ed3ceaSHong Zhang   }
96558ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
96658ed3ceaSHong Zhang 
96740192850SHong Zhang   if (abt->usecoloring) {
96840192850SHong Zhang     MatTransposeColoring matcoloring = abt->matcoloring;
96940192850SHong Zhang     Mat                  Bt_dense,C_dense = abt->ABt_den;
970c8db22f6SHong Zhang 
971b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
97240192850SHong Zhang     Bt_dense = abt->Bt_den;
973b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
974c8db22f6SHong Zhang 
975c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9762128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
977c8db22f6SHong Zhang 
9782128a86cSHong Zhang     /* Recover C from C_dense */
979b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
980c8db22f6SHong Zhang     PetscFunctionReturn(0);
981c8db22f6SHong Zhang   }
982c8db22f6SHong Zhang 
9831fa1209cSHong Zhang   for (i=0; i<cm; i++) {
98481d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
98581d82fe4SHong Zhang     acol = aj + ai[i];
9866973769bSHong Zhang     aval = aa + ai[i];
9871fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9881fa1209cSHong Zhang     ccol = cj + ci[i];
9896973769bSHong Zhang     cval = ca + ci[i];
9901fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
99181d82fe4SHong Zhang       brow = ccol[j];
99281d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
99381d82fe4SHong Zhang       bcol = bj + bi[brow];
9946973769bSHong Zhang       bval = ba + bi[brow];
9956973769bSHong Zhang 
99681d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
99781d82fe4SHong Zhang       nexta = 0; nextb = 0;
99881d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
9997b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
100081d82fe4SHong Zhang         if (nexta == anzi) break;
10017b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
100281d82fe4SHong Zhang         if (nextb == bnzj) break;
100381d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
10046973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
100581d82fe4SHong Zhang           nexta++; nextb++;
100681d82fe4SHong Zhang           flops += 2;
10071fa1209cSHong Zhang         }
10081fa1209cSHong Zhang       }
100981d82fe4SHong Zhang     }
101081d82fe4SHong Zhang   }
101181d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101281d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101381d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10145df89d91SHong Zhang   PetscFunctionReturn(0);
10155df89d91SHong Zhang }
10165df89d91SHong Zhang 
10175df89d91SHong Zhang #undef __FUNCT__
101875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
10190adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10200adebc6cSBarry Smith {
10215df89d91SHong Zhang   PetscErrorCode ierr;
10225df89d91SHong Zhang 
10235df89d91SHong Zhang   PetscFunctionBegin;
10245df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
102507706670SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
102675648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
102707706670SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
10285df89d91SHong Zhang   }
102907706670SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
103075648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
103107706670SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
10325df89d91SHong Zhang   PetscFunctionReturn(0);
10335df89d91SHong Zhang }
10345df89d91SHong Zhang 
10355df89d91SHong Zhang #undef __FUNCT__
103675648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
103775648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10385df89d91SHong Zhang {
1039bc011b1eSHong Zhang   PetscErrorCode ierr;
1040bc011b1eSHong Zhang   Mat            At;
104138baddfdSBarry Smith   PetscInt       *ati,*atj;
1042bc011b1eSHong Zhang 
1043bc011b1eSHong Zhang   PetscFunctionBegin;
1044bc011b1eSHong Zhang   /* create symbolic At */
1045bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10460298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
104733d57670SJed Brown   ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1048bc011b1eSHong Zhang 
1049bc011b1eSHong Zhang   /* get symbolic C=At*B */
1050bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1051bc011b1eSHong Zhang 
1052bc011b1eSHong Zhang   /* clean up */
10536bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1054bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1055bc011b1eSHong Zhang   PetscFunctionReturn(0);
1056bc011b1eSHong Zhang }
1057bc011b1eSHong Zhang 
1058bc011b1eSHong Zhang #undef __FUNCT__
105975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
106075648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1061bc011b1eSHong Zhang {
1062bc011b1eSHong Zhang   PetscErrorCode ierr;
10630fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1064d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1065d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1066d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1067aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1068bc011b1eSHong Zhang 
1069bc011b1eSHong Zhang   PetscFunctionBegin;
1070aa1aec99SHong Zhang   if (!c->a) {
1071854ce69bSBarry Smith     ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
10722205254eSKarl Rupp 
1073aa1aec99SHong Zhang     c->a      = ca;
1074aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1075aa1aec99SHong Zhang   } else {
1076aa1aec99SHong Zhang     ca = c->a;
1077aa1aec99SHong Zhang   }
1078bc011b1eSHong Zhang   /* clear old values in C */
1079bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1080bc011b1eSHong Zhang 
1081bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1082bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1083bc011b1eSHong Zhang     bj   = b->j + bi[i];
1084bc011b1eSHong Zhang     ba   = b->a + bi[i];
1085bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1086bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1087bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1088bc011b1eSHong Zhang       nextb = 0;
10890fbc74f4SHong Zhang       crow  = *aj++;
1090bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1091bc011b1eSHong Zhang       caj   = ca + ci[crow];
1092bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1093bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
10940fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
10950fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1096bc011b1eSHong Zhang           nextb++;
1097bc011b1eSHong Zhang         }
1098bc011b1eSHong Zhang       }
1099bc011b1eSHong Zhang       flops += 2*bnzi;
11000fbc74f4SHong Zhang       aa++;
1101bc011b1eSHong Zhang     }
1102bc011b1eSHong Zhang   }
1103bc011b1eSHong Zhang 
1104bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1105bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1106bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1107bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1108bc011b1eSHong Zhang   PetscFunctionReturn(0);
1109bc011b1eSHong Zhang }
11109123193aSHong Zhang 
11119123193aSHong Zhang #undef __FUNCT__
11129123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
1113150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11149123193aSHong Zhang {
11159123193aSHong Zhang   PetscErrorCode ierr;
11169123193aSHong Zhang 
11179123193aSHong Zhang   PetscFunctionBegin;
11189123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
11193ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11209123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11213ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11229123193aSHong Zhang   }
11233ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11249123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11254614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11269123193aSHong Zhang   PetscFunctionReturn(0);
11279123193aSHong Zhang }
1128edd81eecSMatthew Knepley 
11299123193aSHong Zhang #undef __FUNCT__
11309123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11319123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11329123193aSHong Zhang {
11339123193aSHong Zhang   PetscErrorCode ierr;
11349123193aSHong Zhang 
11359123193aSHong Zhang   PetscFunctionBegin;
11365a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11372205254eSKarl Rupp 
1138d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11399123193aSHong Zhang   PetscFunctionReturn(0);
11409123193aSHong Zhang }
11419123193aSHong Zhang 
11429123193aSHong Zhang #undef __FUNCT__
11439123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11449123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11459123193aSHong Zhang {
1146f32d5d43SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1147612438f1SStefano Zampini   Mat_SeqDense      *bd = (Mat_SeqDense*)B->data;
11489123193aSHong Zhang   PetscErrorCode    ierr;
1149daf57211SHong Zhang   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1150daf57211SHong Zhang   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1151daf57211SHong Zhang   const PetscInt    *aj;
1152612438f1SStefano Zampini   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1153daf57211SHong Zhang   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
11549123193aSHong Zhang 
11559123193aSHong Zhang   PetscFunctionBegin;
1156f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1157612438f1SStefano 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);
1158e32f2f54SBarry 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);
1159e32f2f54SBarry 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);
1160612438f1SStefano Zampini   b = bd->v;
11618c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1162f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1163730858b9SHong Zhang   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1164f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1165f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1166f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1167f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1168f32d5d43SBarry Smith       aj = a->j + a->i[i];
1169f32d5d43SBarry Smith       aa = a->a + a->i[i];
1170f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1171730858b9SHong Zhang         aatmp = aa[j]; ajtmp = aj[j];
1172730858b9SHong Zhang         r1 += aatmp*b1[ajtmp];
1173730858b9SHong Zhang         r2 += aatmp*b2[ajtmp];
1174730858b9SHong Zhang         r3 += aatmp*b3[ajtmp];
1175730858b9SHong Zhang         r4 += aatmp*b4[ajtmp];
11769123193aSHong Zhang       }
1177730858b9SHong Zhang       c1[i] = r1;
1178730858b9SHong Zhang       c2[i] = r2;
1179730858b9SHong Zhang       c3[i] = r3;
1180730858b9SHong Zhang       c4[i] = r4;
1181f32d5d43SBarry Smith     }
1182730858b9SHong Zhang     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1183730858b9SHong Zhang     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1184f32d5d43SBarry Smith   }
1185f32d5d43SBarry Smith   for (; col<cn; col++) {   /* over extra columns of C */
1186f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1187f32d5d43SBarry Smith       r1 = 0.0;
1188f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1189f32d5d43SBarry Smith       aj = a->j + a->i[i];
1190f32d5d43SBarry Smith       aa = a->a + a->i[i];
1191f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1192730858b9SHong Zhang         r1 += aa[j]*b1[aj[j]];
1193f32d5d43SBarry Smith       }
1194730858b9SHong Zhang       c1[i] = r1;
1195f32d5d43SBarry Smith     }
1196f32d5d43SBarry Smith     b1 += bm;
1197730858b9SHong Zhang     c1 += am;
1198f32d5d43SBarry Smith   }
1199dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
12008c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1201f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1202f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1203f32d5d43SBarry Smith   PetscFunctionReturn(0);
1204f32d5d43SBarry Smith }
1205f32d5d43SBarry Smith 
1206f32d5d43SBarry Smith /*
12074324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1208f32d5d43SBarry Smith */
1209f32d5d43SBarry Smith #undef __FUNCT__
1210f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1211f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1212f32d5d43SBarry Smith {
1213f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
121490f5ea3eSStefano Zampini   Mat_SeqDense   *bd = (Mat_SeqDense*)B->data;
1215f32d5d43SBarry Smith   PetscErrorCode ierr;
1216dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1217dd6ea824SBarry Smith   MatScalar      *aa;
121890f5ea3eSStefano Zampini   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
12194324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1220f32d5d43SBarry Smith 
1221f32d5d43SBarry Smith   PetscFunctionBegin;
1222f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
122390f5ea3eSStefano Zampini   b = bd->v;
12248c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1225f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12264324174eSBarry Smith 
12274324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12284324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12294324174eSBarry Smith       colam = col*am;
12304324174eSBarry Smith       arm   = a->compressedrow.nrows;
12314324174eSBarry Smith       ii    = a->compressedrow.i;
12324324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12334324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12344324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12354324174eSBarry Smith         n  = ii[i+1] - ii[i];
12364324174eSBarry Smith         aj = a->j + ii[i];
12374324174eSBarry Smith         aa = a->a + ii[i];
12384324174eSBarry Smith         for (j=0; j<n; j++) {
12394324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12404324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12414324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12424324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12434324174eSBarry Smith         }
12444324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12454324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12464324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12474324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12484324174eSBarry Smith       }
12494324174eSBarry Smith       b1 += bm4;
12504324174eSBarry Smith       b2 += bm4;
12514324174eSBarry Smith       b3 += bm4;
12524324174eSBarry Smith       b4 += bm4;
12534324174eSBarry Smith     }
12544324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12554324174eSBarry Smith       colam = col*am;
12564324174eSBarry Smith       arm   = a->compressedrow.nrows;
12574324174eSBarry Smith       ii    = a->compressedrow.i;
12584324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12594324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12604324174eSBarry Smith         r1 = 0.0;
12614324174eSBarry Smith         n  = ii[i+1] - ii[i];
12624324174eSBarry Smith         aj = a->j + ii[i];
12634324174eSBarry Smith         aa = a->a + ii[i];
12644324174eSBarry Smith 
12654324174eSBarry Smith         for (j=0; j<n; j++) {
12664324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
12674324174eSBarry Smith         }
1268a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
12694324174eSBarry Smith       }
12704324174eSBarry Smith       b1 += bm;
12714324174eSBarry Smith     }
12724324174eSBarry Smith   } else {
1273f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1274f32d5d43SBarry Smith       colam = col*am;
1275f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1276f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1277f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1278f32d5d43SBarry Smith         aj = a->j + a->i[i];
1279f32d5d43SBarry Smith         aa = a->a + a->i[i];
1280f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1281f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1282f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1283f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1284f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1285f32d5d43SBarry Smith         }
1286f32d5d43SBarry Smith         c[colam + i]       += r1;
1287f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1288f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1289f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1290f32d5d43SBarry Smith       }
1291f32d5d43SBarry Smith       b1 += bm4;
1292f32d5d43SBarry Smith       b2 += bm4;
1293f32d5d43SBarry Smith       b3 += bm4;
1294f32d5d43SBarry Smith       b4 += bm4;
1295f32d5d43SBarry Smith     }
1296f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1297a2ea699eSBarry Smith       colam = col*am;
1298f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1299f32d5d43SBarry Smith         r1 = 0.0;
1300f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1301f32d5d43SBarry Smith         aj = a->j + a->i[i];
1302f32d5d43SBarry Smith         aa = a->a + a->i[i];
1303f32d5d43SBarry Smith 
1304f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1305f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1306f32d5d43SBarry Smith         }
1307a2ea699eSBarry Smith         c[colam + i] += r1;
1308f32d5d43SBarry Smith       }
1309f32d5d43SBarry Smith       b1 += bm;
1310f32d5d43SBarry Smith     }
13114324174eSBarry Smith   }
1312dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13138c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13149123193aSHong Zhang   PetscFunctionReturn(0);
13159123193aSHong Zhang }
1316b1683b59SHong Zhang 
1317b1683b59SHong Zhang #undef __FUNCT__
1318b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1319b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1320c8db22f6SHong Zhang {
1321c8db22f6SHong Zhang   PetscErrorCode ierr;
13222f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13232f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13242f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13252f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13262f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13272f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1328c8db22f6SHong Zhang 
1329c8db22f6SHong Zhang   PetscFunctionBegin;
13302f5f1f90SHong Zhang   btval_den=btdense->v;
13312f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13322f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13332f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13342f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1335d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13362f5f1f90SHong Zhang       btcol = bj + bi[col];
13372f5f1f90SHong Zhang       btval = ba + bi[col];
13382f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1339d2853540SHong Zhang       for (j=0; j<anz; j++) {
13402f5f1f90SHong Zhang         brow            = btcol[j];
13412f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1342c8db22f6SHong Zhang       }
1343c8db22f6SHong Zhang     }
13442f5f1f90SHong Zhang     btval_den += m;
1345c8db22f6SHong Zhang   }
1346c8db22f6SHong Zhang   PetscFunctionReturn(0);
1347c8db22f6SHong Zhang }
1348c8db22f6SHong Zhang 
1349c8db22f6SHong Zhang #undef __FUNCT__
1350b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1351b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13528972f759SHong Zhang {
13538972f759SHong Zhang   PetscErrorCode ierr;
1354b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
1355077f23c2SHong Zhang   PetscScalar    *ca_den,*ca_den_ptr,*ca=csp->a;
1356f99a636bSHong Zhang   PetscInt       k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1357e88777f2SHong Zhang   PetscInt       brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1358077f23c2SHong Zhang   PetscInt       nrows,*row,*idx;
1359077f23c2SHong Zhang   PetscInt       *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
13608972f759SHong Zhang 
13618972f759SHong Zhang   PetscFunctionBegin;
1362a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1363f99a636bSHong Zhang 
1364077f23c2SHong Zhang   if (brows > 0) {
1365077f23c2SHong Zhang     PetscInt *lstart,row_end,row_start;
1366980ae229SHong Zhang     lstart = matcoloring->lstart;
1367eeb4fd02SHong Zhang     ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr);
1368eeb4fd02SHong Zhang 
1369077f23c2SHong Zhang     row_end = brows;
1370eeb4fd02SHong Zhang     if (row_end > m) row_end = m;
1371077f23c2SHong Zhang     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1372077f23c2SHong Zhang       ca_den_ptr = ca_den;
1373980ae229SHong Zhang       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1374eeb4fd02SHong Zhang         nrows = matcoloring->nrows[k];
1375eeb4fd02SHong Zhang         row   = rows  + colorforrow[k];
1376eeb4fd02SHong Zhang         idx   = den2sp + colorforrow[k];
1377eeb4fd02SHong Zhang         for (l=lstart[k]; l<nrows; l++) {
1378eeb4fd02SHong Zhang           if (row[l] >= row_end) {
1379eeb4fd02SHong Zhang             lstart[k] = l;
1380eeb4fd02SHong Zhang             break;
1381eeb4fd02SHong Zhang           } else {
1382077f23c2SHong Zhang             ca[idx[l]] = ca_den_ptr[row[l]];
1383eeb4fd02SHong Zhang           }
1384eeb4fd02SHong Zhang         }
1385077f23c2SHong Zhang         ca_den_ptr += m;
1386eeb4fd02SHong Zhang       }
1387077f23c2SHong Zhang       row_end += brows;
1388eeb4fd02SHong Zhang       if (row_end > m) row_end = m;
1389eeb4fd02SHong Zhang     }
1390077f23c2SHong Zhang   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1391077f23c2SHong Zhang     ca_den_ptr = ca_den;
1392077f23c2SHong Zhang     for (k=0; k<ncolors; k++) {
1393077f23c2SHong Zhang       nrows = matcoloring->nrows[k];
1394077f23c2SHong Zhang       row   = rows  + colorforrow[k];
1395077f23c2SHong Zhang       idx   = den2sp + colorforrow[k];
1396077f23c2SHong Zhang       for (l=0; l<nrows; l++) {
1397077f23c2SHong Zhang         ca[idx[l]] = ca_den_ptr[row[l]];
1398077f23c2SHong Zhang       }
1399077f23c2SHong Zhang       ca_den_ptr += m;
1400077f23c2SHong Zhang     }
1401f99a636bSHong Zhang   }
1402f99a636bSHong Zhang 
1403a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
1404f99a636bSHong Zhang #if defined(PETSC_USE_INFO)
1405077f23c2SHong Zhang   if (matcoloring->brows > 0) {
1406f99a636bSHong Zhang     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1407e88777f2SHong Zhang   } else {
1408077f23c2SHong Zhang     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1409e88777f2SHong Zhang   }
1410f99a636bSHong Zhang #endif
14118972f759SHong Zhang   PetscFunctionReturn(0);
14128972f759SHong Zhang }
14138972f759SHong Zhang 
14148972f759SHong Zhang #undef __FUNCT__
1415b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1416b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1417b1683b59SHong Zhang {
1418b1683b59SHong Zhang   PetscErrorCode ierr;
1419e88777f2SHong Zhang   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
14201a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1421b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1422b1683b59SHong Zhang   IS             *isa;
1423b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1424e88777f2SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1425e88777f2SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i,brows;
1426e88777f2SHong Zhang   PetscBool      flg;
1427b1683b59SHong Zhang 
1428b1683b59SHong Zhang   PetscFunctionBegin;
1429b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1430e99be685SHong Zhang 
14317ecccc15SHong Zhang   /* bs >1 is not being tested yet! */
1432e88777f2SHong Zhang   Nbs       = mat->cmap->N/bs;
1433b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1434e88777f2SHong Zhang   c->N      = Nbs;
1435e88777f2SHong Zhang   c->m      = c->M;
1436b1683b59SHong Zhang   c->rstart = 0;
1437e88777f2SHong Zhang   c->brows  = 100;
1438b1683b59SHong Zhang 
1439b1683b59SHong Zhang   c->ncolors = nis;
1440dcca6d9dSJed Brown   ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr);
1441854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr);
1442854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr);
1443e88777f2SHong Zhang 
1444e88777f2SHong Zhang   brows = c->brows;
1445c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1446e88777f2SHong Zhang   if (flg) c->brows = brows;
1447eeb4fd02SHong Zhang   if (brows > 0) {
1448854ce69bSBarry Smith     ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr);
1449e88777f2SHong Zhang   }
14502205254eSKarl Rupp 
1451d2853540SHong Zhang   colorforrow[0] = 0;
1452d2853540SHong Zhang   rows_i         = rows;
1453f99a636bSHong Zhang   den2sp_i       = den2sp;
1454b1683b59SHong Zhang 
1455854ce69bSBarry Smith   ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr);
1456854ce69bSBarry Smith   ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr);
14572205254eSKarl Rupp 
1458d2853540SHong Zhang   colorforcol[0] = 0;
1459d2853540SHong Zhang   columns_i      = columns;
1460d2853540SHong Zhang 
1461077f23c2SHong Zhang   /* get column-wise storage of mat */
1462077f23c2SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1463b1683b59SHong Zhang 
1464b28d80bdSHong Zhang   cm   = c->m;
1465854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr);
1466854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr);
1467f99a636bSHong Zhang   for (i=0; i<nis; i++) { /* loop over color */
1468b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1469b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
14702205254eSKarl Rupp 
1471b1683b59SHong Zhang     c->ncolumns[i] = n;
1472b1683b59SHong Zhang     if (n) {
1473d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1474b1683b59SHong Zhang     }
1475d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1476d2853540SHong Zhang     columns_i       += n;
1477b1683b59SHong Zhang 
1478b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1479e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1480e99be685SHong Zhang 
1481b7caf3d6SHong Zhang     for (j=0; j<n; j++) { /* loop over columns*/
1482b1683b59SHong Zhang       col     = is[j];
1483d2853540SHong Zhang       row_idx = cj + ci[col];
1484b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1485b7caf3d6SHong Zhang       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1486e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1487d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1488b1683b59SHong Zhang       }
1489b1683b59SHong Zhang     }
1490b1683b59SHong Zhang     /* count the number of hits */
1491b1683b59SHong Zhang     nrows = 0;
1492e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1493b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1494b1683b59SHong Zhang     }
1495b1683b59SHong Zhang     c->nrows[i]      = nrows;
1496d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1497d2853540SHong Zhang 
1498b1683b59SHong Zhang     nrows = 0;
1499b7caf3d6SHong Zhang     for (j=0; j<cm; j++) { /* loop over rows */
1500b1683b59SHong Zhang       if (rowhit[j]) {
1501d2853540SHong Zhang         rows_i[nrows]   = j;
150212b89a43SHong Zhang         den2sp_i[nrows] = idxhit[j];
1503b1683b59SHong Zhang         nrows++;
1504b1683b59SHong Zhang       }
1505b1683b59SHong Zhang     }
1506e88777f2SHong Zhang     den2sp_i += nrows;
1507077f23c2SHong Zhang 
1508b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1509f99a636bSHong Zhang     rows_i += nrows;
1510b1683b59SHong Zhang   }
15110298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1512b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1513b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1514d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1515b28d80bdSHong Zhang 
1516d2853540SHong Zhang   c->colorforrow = colorforrow;
1517d2853540SHong Zhang   c->rows        = rows;
1518f99a636bSHong Zhang   c->den2sp      = den2sp;
1519d2853540SHong Zhang   c->colorforcol = colorforcol;
1520d2853540SHong Zhang   c->columns     = columns;
15212205254eSKarl Rupp 
1522f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1523b1683b59SHong Zhang   PetscFunctionReturn(0);
1524b1683b59SHong Zhang }
1525