xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 4d478ae76771a501e2c3958c67bfeb93cdb5fe87)
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>
1107475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h>
1226be7272SHong Zhang #include <petsctime.h>
137bab7c10SHong Zhang 
1458cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
15cd093f1dSJed Brown 
166284ec50SHong Zhang #undef __FUNCT__
175c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1838baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1938baddfdSBarry Smith {
20dfbe8321SBarry Smith   PetscErrorCode ierr;
2158cf0668SJed Brown   PetscBool      scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE,btheap = PETSC_FALSE,llcondensed = PETSC_FALSE;
225c66b693SKris Buschelman 
235c66b693SKris Buschelman   PetscFunctionBegin;
2426be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
25152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
260298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,NULL);CHKERRQ(ierr);
270298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,NULL);CHKERRQ(ierr);
280298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,NULL);CHKERRQ(ierr);
290298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_btheap","Use btheap implementation of symbolic factorization C=A*B","",btheap,&btheap,NULL);CHKERRQ(ierr);
3058cf0668SJed Brown     ierr = PetscOptionsBool("-matmatmult_llcondensed","Use LLCondensed to for symbolic C=A*B","",llcondensed,&llcondensed,NULL);CHKERRQ(ierr);
31d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
323ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
333c50cad2SHong Zhang     if (scalable_fast) {
343c50cad2SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
35aacf854cSBarry Smith     } else if (scalable) {
36aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
370ced3a2bSJed Brown     } else if (heap) {
380ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
398a07c6f1SJed Brown     } else if (btheap) {
408a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
4158cf0668SJed Brown     } else if (llcondensed) {
4258cf0668SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
4325023028SHong Zhang     } else {
4426be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
4525023028SHong Zhang     }
463ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
4726be0446SHong Zhang   }
485c913ed7SHong Zhang 
493ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
5001e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
513ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
525c66b693SKris Buschelman   PetscFunctionReturn(0);
535c66b693SKris Buschelman }
541c24bd37SHong Zhang 
55e9e4536cSHong Zhang #undef __FUNCT__
5658cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed"
5758cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
58b561aa9dSHong Zhang {
59b561aa9dSHong Zhang   PetscErrorCode     ierr;
60b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
61971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
62c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
63b561aa9dSHong Zhang   PetscReal          afill;
64fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
65fb908031SHong Zhang   PetscBT            lnkbt;
660298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
67b561aa9dSHong Zhang 
68b561aa9dSHong Zhang   PetscFunctionBegin;
69bd958071SHong Zhang   /* Get ci and cj */
70bd958071SHong Zhang   /*---------------*/
71fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
72fb908031SHong Zhang   /* free space for accumulating nonzero column info */
73fb908031SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
74fb908031SHong Zhang   ci[0] = 0;
75fb908031SHong Zhang 
76fb908031SHong Zhang   /* create and initialize a linked list */
77fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
78fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
79fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
80fb908031SHong Zhang 
81fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
82fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
832205254eSKarl Rupp 
84fb908031SHong Zhang   current_space = free_space;
85fb908031SHong Zhang 
86fb908031SHong Zhang   /* Determine ci and cj */
87fb908031SHong Zhang   for (i=0; i<am; i++) {
88fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
89fb908031SHong Zhang     aj   = a->j + ai[i];
90fb908031SHong Zhang     for (j=0; j<anzi; j++) {
91fb908031SHong Zhang       brow = aj[j];
92fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
93fb908031SHong Zhang       bj   = b->j + bi[brow];
94fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
95fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
96fb908031SHong Zhang     }
97fb908031SHong Zhang     cnzi = lnk[0];
98fb908031SHong Zhang 
99fb908031SHong Zhang     /* If free space is not available, make more free space */
100fb908031SHong Zhang     /* Double the amount of total space in the list */
101fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
102fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
103fb908031SHong Zhang       ndouble++;
104fb908031SHong Zhang     }
105fb908031SHong Zhang 
106fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
107fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1082205254eSKarl Rupp 
109fb908031SHong Zhang     current_space->array           += cnzi;
110fb908031SHong Zhang     current_space->local_used      += cnzi;
111fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1122205254eSKarl Rupp 
113fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
114fb908031SHong Zhang   }
115fb908031SHong Zhang 
116fb908031SHong Zhang   /* Column indices are in the list of free space */
117fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
118fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
119fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
120fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
121fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
122b561aa9dSHong Zhang 
12326be0446SHong Zhang   /* put together the new symbolic matrix */
124ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1252205254eSKarl Rupp 
126a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
127a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
12858c24d83SHong Zhang 
12958c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
13058c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
13158c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
132aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
133e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
13458c24d83SHong Zhang   c->nonew                  = 0;
135002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1360b7e3e3dSHong Zhang 
1377212da7cSHong Zhang   /* set MatInfo */
1387212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
139f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1407212da7cSHong Zhang   c->maxnz                     = ci[am];
1417212da7cSHong Zhang   c->nz                        = ci[am];
142fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1437212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1447212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1457212da7cSHong Zhang 
1467212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1477212da7cSHong Zhang   if (ci[am]) {
148fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
149f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
150f2b054eeSHong Zhang   } else {
151f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
152be0fcf8dSHong Zhang   }
153f2b054eeSHong Zhang #endif
15458c24d83SHong Zhang   PetscFunctionReturn(0);
15558c24d83SHong Zhang }
156d50806bdSBarry Smith 
15726be0446SHong Zhang #undef __FUNCT__
15826be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
159dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
160d50806bdSBarry Smith {
161dfbe8321SBarry Smith   PetscErrorCode ierr;
162d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
163d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
164d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
165d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
16638baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
167c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
168519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
169aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
170aa1aec99SHong Zhang   PetscScalar    *ab_dense;
171d50806bdSBarry Smith 
172d50806bdSBarry Smith   PetscFunctionBegin;
173aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
174aa1aec99SHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
175aa1aec99SHong Zhang     c->a      = ca;
176aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
177aa1aec99SHong Zhang 
178aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
179aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
1802205254eSKarl Rupp 
181aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
182aa1aec99SHong Zhang   } else {
183aa1aec99SHong Zhang     ca       = c->a;
184aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
185aa1aec99SHong Zhang   }
186aa1aec99SHong Zhang 
187c124e916SHong Zhang   /* clean old values in C */
188c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
189d50806bdSBarry Smith   /* Traverse A row-wise. */
190d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
191d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
192d50806bdSBarry Smith   for (i=0; i<am; i++) {
193d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
194d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
195519eb980SHong Zhang       brow = aj[j];
196d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
197d50806bdSBarry Smith       bjj  = bj + bi[brow];
198d50806bdSBarry Smith       baj  = ba + bi[brow];
199519eb980SHong Zhang       /* perform dense axpy */
20036ec6d2dSHong Zhang       valtmp = aa[j];
201519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
20236ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
203519eb980SHong Zhang       }
204519eb980SHong Zhang       flops += 2*bnzi;
205519eb980SHong Zhang     }
206c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
207c58ca2e3SHong Zhang 
208c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
209519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
210519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
211519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
212519eb980SHong Zhang     }
213519eb980SHong Zhang     flops += cnzi;
214519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
215519eb980SHong Zhang   }
216c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
217c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
218c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
219c58ca2e3SHong Zhang   PetscFunctionReturn(0);
220c58ca2e3SHong Zhang }
221c58ca2e3SHong Zhang 
222c58ca2e3SHong Zhang #undef __FUNCT__
22325023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
22425023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
225c58ca2e3SHong Zhang {
226c58ca2e3SHong Zhang   PetscErrorCode ierr;
227c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
228c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
229c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
230c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
231c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
232c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
233c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
23436ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
235c58ca2e3SHong Zhang   PetscInt       nextb;
236c58ca2e3SHong Zhang 
237c58ca2e3SHong Zhang   PetscFunctionBegin;
238c58ca2e3SHong Zhang   /* clean old values in C */
239c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
240c58ca2e3SHong Zhang   /* Traverse A row-wise. */
241c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
242c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
243519eb980SHong Zhang   for (i=0; i<am; i++) {
244519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
245519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
246519eb980SHong Zhang     for (j=0; j<anzi; j++) {
247519eb980SHong Zhang       brow = aj[j];
248519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
249519eb980SHong Zhang       bjj  = bj + bi[brow];
250519eb980SHong Zhang       baj  = ba + bi[brow];
251519eb980SHong Zhang       /* perform sparse axpy */
25236ec6d2dSHong Zhang       valtmp = aa[j];
253c124e916SHong Zhang       nextb  = 0;
254c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
255c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
25636ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
257c124e916SHong Zhang         }
258d50806bdSBarry Smith       }
259d50806bdSBarry Smith       flops += 2*bnzi;
260d50806bdSBarry Smith     }
261519eb980SHong Zhang     aj += anzi; aa += anzi;
262519eb980SHong Zhang     cj += cnzi; ca += cnzi;
263d50806bdSBarry Smith   }
264c58ca2e3SHong Zhang 
265716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
266716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
267d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
268d50806bdSBarry Smith   PetscFunctionReturn(0);
269d50806bdSBarry Smith }
270bc011b1eSHong Zhang 
271e9e4536cSHong Zhang #undef __FUNCT__
2723c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2733c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
27425296bd5SBarry Smith {
27525296bd5SBarry Smith   PetscErrorCode     ierr;
27625296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
27725296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2783c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
27925296bd5SBarry Smith   MatScalar          *ca;
28025296bd5SBarry Smith   PetscReal          afill;
2813c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2820298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
28325296bd5SBarry Smith 
28425296bd5SBarry Smith   PetscFunctionBegin;
2853c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2863c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2873c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2883c50cad2SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2893c50cad2SHong Zhang   ci[0] = 0;
2903c50cad2SHong Zhang 
2913c50cad2SHong Zhang   /* create and initialize a linked list */
2923c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2933c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2943c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
2953c50cad2SHong Zhang 
2963c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
2973c50cad2SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
2983c50cad2SHong Zhang   current_space = free_space;
2993c50cad2SHong Zhang 
3003c50cad2SHong Zhang   /* Determine ci and cj */
3013c50cad2SHong Zhang   for (i=0; i<am; i++) {
3023c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3033c50cad2SHong Zhang     aj   = a->j + ai[i];
3043c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3053c50cad2SHong Zhang       brow = aj[j];
3063c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3073c50cad2SHong Zhang       bj   = b->j + bi[brow];
3083c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3093c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3103c50cad2SHong Zhang     }
3113c50cad2SHong Zhang     cnzi = lnk[1];
3123c50cad2SHong Zhang 
3133c50cad2SHong Zhang     /* If free space is not available, make more free space */
3143c50cad2SHong Zhang     /* Double the amount of total space in the list */
3153c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3163c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3173c50cad2SHong Zhang       ndouble++;
3183c50cad2SHong Zhang     }
3193c50cad2SHong Zhang 
3203c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3213c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3222205254eSKarl Rupp 
3233c50cad2SHong Zhang     current_space->array           += cnzi;
3243c50cad2SHong Zhang     current_space->local_used      += cnzi;
3253c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3262205254eSKarl Rupp 
3273c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3283c50cad2SHong Zhang   }
3293c50cad2SHong Zhang 
3303c50cad2SHong Zhang   /* Column indices are in the list of free space */
3313c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3323c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3333c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3343c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3353c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
33625296bd5SBarry Smith 
33725296bd5SBarry Smith   /* Allocate space for ca */
33825296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
33925296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
34025296bd5SBarry Smith 
34125296bd5SBarry Smith   /* put together the new symbolic matrix */
342ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
3432205254eSKarl Rupp 
344a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
345a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
34625296bd5SBarry Smith 
34725296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
34825296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
34925296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
35025296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
35125296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
35225296bd5SBarry Smith   c->nonew   = 0;
3532205254eSKarl Rupp 
35425296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
35525296bd5SBarry Smith 
35625296bd5SBarry Smith   /* set MatInfo */
35725296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
35825296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
35925296bd5SBarry Smith   c->maxnz                     = ci[am];
36025296bd5SBarry Smith   c->nz                        = ci[am];
3613c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
36225296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
36325296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
36425296bd5SBarry Smith 
36525296bd5SBarry Smith #if defined(PETSC_USE_INFO)
36625296bd5SBarry Smith   if (ci[am]) {
3673c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
36825296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
36925296bd5SBarry Smith   } else {
37025296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
37125296bd5SBarry Smith   }
37225296bd5SBarry Smith #endif
37325296bd5SBarry Smith   PetscFunctionReturn(0);
37425296bd5SBarry Smith }
37525296bd5SBarry Smith 
37625296bd5SBarry Smith 
37725296bd5SBarry Smith #undef __FUNCT__
37825023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
37925023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
380e9e4536cSHong Zhang {
381e9e4536cSHong Zhang   PetscErrorCode     ierr;
382e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
383bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
38425c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
385e9e4536cSHong Zhang   MatScalar          *ca;
386e9e4536cSHong Zhang   PetscReal          afill;
387bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
3880298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
389e9e4536cSHong Zhang 
390e9e4536cSHong Zhang   PetscFunctionBegin;
391bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
392bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
393bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
394bd958071SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
395bd958071SHong Zhang   ci[0] = 0;
396bd958071SHong Zhang 
397bd958071SHong Zhang   /* create and initialize a linked list */
398bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
399bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
400bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
401bd958071SHong Zhang 
402bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
403bd958071SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
404bd958071SHong Zhang   current_space = free_space;
405bd958071SHong Zhang 
406bd958071SHong Zhang   /* Determine ci and cj */
407bd958071SHong Zhang   for (i=0; i<am; i++) {
408bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
409bd958071SHong Zhang     aj   = a->j + ai[i];
410bd958071SHong Zhang     for (j=0; j<anzi; j++) {
411bd958071SHong Zhang       brow = aj[j];
412bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
413bd958071SHong Zhang       bj   = b->j + bi[brow];
414bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
415bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
416bd958071SHong Zhang     }
417bd958071SHong Zhang     cnzi = lnk[0];
418bd958071SHong Zhang 
419bd958071SHong Zhang     /* If free space is not available, make more free space */
420bd958071SHong Zhang     /* Double the amount of total space in the list */
421bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
422bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
423bd958071SHong Zhang       ndouble++;
424bd958071SHong Zhang     }
425bd958071SHong Zhang 
426bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
427bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4282205254eSKarl Rupp 
429bd958071SHong Zhang     current_space->array           += cnzi;
430bd958071SHong Zhang     current_space->local_used      += cnzi;
431bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4322205254eSKarl Rupp 
433bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
434bd958071SHong Zhang   }
435bd958071SHong Zhang 
436bd958071SHong Zhang   /* Column indices are in the list of free space */
437bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
438bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
439bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
440bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
441bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
442e9e4536cSHong Zhang 
443e9e4536cSHong Zhang   /* Allocate space for ca */
444bd958071SHong Zhang   /*-----------------------*/
445e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
446e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
447e9e4536cSHong Zhang 
448e9e4536cSHong Zhang   /* put together the new symbolic matrix */
449ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
4502205254eSKarl Rupp 
451a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
452a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
453e9e4536cSHong Zhang 
454e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
455e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
456e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
457e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
458e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
459e9e4536cSHong Zhang   c->nonew   = 0;
4602205254eSKarl Rupp 
46125023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
462e9e4536cSHong Zhang 
463e9e4536cSHong Zhang   /* set MatInfo */
464e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
465e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
466e9e4536cSHong Zhang   c->maxnz                     = ci[am];
467e9e4536cSHong Zhang   c->nz                        = ci[am];
468bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
469e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
470e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
471e9e4536cSHong Zhang 
472e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
473e9e4536cSHong Zhang   if (ci[am]) {
474bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
475e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
476e9e4536cSHong Zhang   } else {
477e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
478e9e4536cSHong Zhang   }
479e9e4536cSHong Zhang #endif
480e9e4536cSHong Zhang   PetscFunctionReturn(0);
481e9e4536cSHong Zhang }
482e9e4536cSHong Zhang 
4830ced3a2bSJed Brown #undef __FUNCT__
4840ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
4850ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
4860ced3a2bSJed Brown {
4870ced3a2bSJed Brown   PetscErrorCode     ierr;
4880ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4890ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
4900ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
4910ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
4920ced3a2bSJed Brown   PetscReal          afill;
4930ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
4940298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
4950ced3a2bSJed Brown   PetscHeap          h;
4960ced3a2bSJed Brown 
4970ced3a2bSJed Brown   PetscFunctionBegin;
498cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
4990ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
5000ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
5010ced3a2bSJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5020ced3a2bSJed Brown   ci[0] = 0;
5030ced3a2bSJed Brown 
5040ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
5050ced3a2bSJed Brown   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5060ced3a2bSJed Brown   current_space = free_space;
5070ced3a2bSJed Brown 
5080ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
5090ced3a2bSJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
5100ced3a2bSJed Brown 
5110ced3a2bSJed Brown   /* Determine ci and cj */
5120ced3a2bSJed Brown   for (i=0; i<am; i++) {
5130ced3a2bSJed 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 */
5140ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5150ced3a2bSJed Brown     ci[i+1] = ci[i];
5160ced3a2bSJed Brown     /* Populate the min heap */
5170ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5180ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5190ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5200ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5210ced3a2bSJed Brown       }
5220ced3a2bSJed Brown     }
5230ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5240ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5250ced3a2bSJed Brown     while (j >= 0) {
5260ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
5270ced3a2bSJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
5280ced3a2bSJed Brown         ndouble++;
5290ced3a2bSJed Brown       }
5300ced3a2bSJed Brown       *(current_space->array++) = col;
5310ced3a2bSJed Brown       current_space->local_used++;
5320ced3a2bSJed Brown       current_space->local_remaining--;
5330ced3a2bSJed Brown       ci[i+1]++;
5340ced3a2bSJed Brown 
5350ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5360ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5370ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5380ced3a2bSJed Brown         PetscInt j2,col2;
5390ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5400ced3a2bSJed Brown         if (col2 != col) break;
5410ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5420ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5430ced3a2bSJed Brown       }
5440ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5450ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5460ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5470ced3a2bSJed Brown     }
5480ced3a2bSJed Brown   }
5490ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5500ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5510ced3a2bSJed Brown 
5520ced3a2bSJed Brown   /* Column indices are in the list of free space */
5530ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5540ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
5550ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5560ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5570ced3a2bSJed Brown 
5580ced3a2bSJed Brown   /* put together the new symbolic matrix */
559ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
5602205254eSKarl Rupp 
5610ced3a2bSJed Brown   (*C)->rmap->bs = A->rmap->bs;
5620ced3a2bSJed Brown   (*C)->cmap->bs = B->cmap->bs;
5630ced3a2bSJed Brown 
5640ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5650ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5660ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5670ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5680ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5690ced3a2bSJed Brown   c->nonew   = 0;
57026fbe8dcSKarl Rupp 
57189d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5720ced3a2bSJed Brown 
5730ced3a2bSJed Brown   /* set MatInfo */
5740ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5750ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5760ced3a2bSJed Brown   c->maxnz                     = ci[am];
5770ced3a2bSJed Brown   c->nz                        = ci[am];
5780ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5790ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5800ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5810ced3a2bSJed Brown 
5820ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5830ced3a2bSJed Brown   if (ci[am]) {
5840ced3a2bSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
5850ced3a2bSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
5860ced3a2bSJed Brown   } else {
5870ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
5880ced3a2bSJed Brown   }
5890ced3a2bSJed Brown #endif
5900ced3a2bSJed Brown   PetscFunctionReturn(0);
5910ced3a2bSJed Brown }
592e9e4536cSHong Zhang 
5938a07c6f1SJed Brown #undef __FUNCT__
5948a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
5958a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
5968a07c6f1SJed Brown {
5978a07c6f1SJed Brown   PetscErrorCode     ierr;
5988a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5998a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
6008a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6018a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6028a07c6f1SJed Brown   PetscReal          afill;
6038a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6040298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6058a07c6f1SJed Brown   PetscHeap          h;
6068a07c6f1SJed Brown   PetscBT            bt;
6078a07c6f1SJed Brown 
6088a07c6f1SJed Brown   PetscFunctionBegin;
609cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6108a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6118a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
6128a07c6f1SJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6138a07c6f1SJed Brown   ci[0] = 0;
6148a07c6f1SJed Brown 
6158a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
6168a07c6f1SJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6172205254eSKarl Rupp 
6188a07c6f1SJed Brown   current_space = free_space;
6198a07c6f1SJed Brown 
6208a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
6218a07c6f1SJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
6228a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6238a07c6f1SJed Brown 
6248a07c6f1SJed Brown   /* Determine ci and cj */
6258a07c6f1SJed Brown   for (i=0; i<am; i++) {
6268a07c6f1SJed 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 */
6278a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6288a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6298a07c6f1SJed Brown     ci[i+1] = ci[i];
6308a07c6f1SJed Brown     /* Populate the min heap */
6318a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6328a07c6f1SJed Brown       PetscInt brow = acol[j];
6338a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6348a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6358a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6368a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6378a07c6f1SJed Brown           bb[j]++;
6388a07c6f1SJed Brown           break;
6398a07c6f1SJed Brown         }
6408a07c6f1SJed Brown       }
6418a07c6f1SJed Brown     }
6428a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6438a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6448a07c6f1SJed Brown     while (j >= 0) {
6458a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6460298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
6478a07c6f1SJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
6488a07c6f1SJed Brown         ndouble++;
6498a07c6f1SJed Brown       }
6508a07c6f1SJed Brown       *(current_space->array++) = col;
6518a07c6f1SJed Brown       current_space->local_used++;
6528a07c6f1SJed Brown       current_space->local_remaining--;
6538a07c6f1SJed Brown       ci[i+1]++;
6548a07c6f1SJed Brown 
6558a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6568a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6578a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6588a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6598a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6608a07c6f1SJed Brown           bb[j]++;
6618a07c6f1SJed Brown           break;
6628a07c6f1SJed Brown         }
6638a07c6f1SJed Brown       }
6648a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6658a07c6f1SJed Brown     }
6668a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6678a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6688a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6698a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6708a07c6f1SJed Brown     }
6718a07c6f1SJed Brown   }
6728a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6738a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6748a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6758a07c6f1SJed Brown 
6768a07c6f1SJed Brown   /* Column indices are in the list of free space */
6778a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6788a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
6798a07c6f1SJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6808a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6818a07c6f1SJed Brown 
6828a07c6f1SJed Brown   /* put together the new symbolic matrix */
683ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
6842205254eSKarl Rupp 
6858a07c6f1SJed Brown   (*C)->rmap->bs = A->rmap->bs;
6868a07c6f1SJed Brown   (*C)->cmap->bs = B->cmap->bs;
6878a07c6f1SJed Brown 
6888a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6898a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6908a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
6918a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
6928a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
6938a07c6f1SJed Brown   c->nonew   = 0;
69426fbe8dcSKarl Rupp 
69589d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
6968a07c6f1SJed Brown 
6978a07c6f1SJed Brown   /* set MatInfo */
6988a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6998a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
7008a07c6f1SJed Brown   c->maxnz                     = ci[am];
7018a07c6f1SJed Brown   c->nz                        = ci[am];
7028a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7038a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7048a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7058a07c6f1SJed Brown 
7068a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7078a07c6f1SJed Brown   if (ci[am]) {
7088a07c6f1SJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
7098a07c6f1SJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
7108a07c6f1SJed Brown   } else {
7118a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7128a07c6f1SJed Brown   }
7138a07c6f1SJed Brown #endif
7148a07c6f1SJed Brown   PetscFunctionReturn(0);
7158a07c6f1SJed Brown }
7168a07c6f1SJed Brown 
717cd093f1dSJed Brown #undef __FUNCT__
71858cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
719cd093f1dSJed Brown /* concatenate unique entries and then sort */
72058cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
721cd093f1dSJed Brown {
722cd093f1dSJed Brown   PetscErrorCode     ierr;
723cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
724cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
725cd093f1dSJed Brown   PetscInt           *ci,*cj;
726cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
727cd093f1dSJed Brown   PetscReal          afill;
728cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
729cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
730cd093f1dSJed Brown   char               *seen;
731cd093f1dSJed Brown 
732cd093f1dSJed Brown   PetscFunctionBegin;
733cd093f1dSJed Brown   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
734cd093f1dSJed Brown   ci[0] = 0;
735cd093f1dSJed Brown 
736cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
737cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
738cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
739cd093f1dSJed Brown   ierr = PetscMalloc(bn*sizeof(char),&seen);CHKERRQ(ierr);
740cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
741cd093f1dSJed Brown 
742cd093f1dSJed Brown   /* Determine ci and cj */
743cd093f1dSJed Brown   for (i=0; i<am; i++) {
744cd093f1dSJed 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 */
745cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
746cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
747cd093f1dSJed Brown     /* Pack segrow */
748cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
749cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
750cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
751cd093f1dSJed Brown         PetscInt bcol = bj[k];
752cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
753cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
754cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
755cd093f1dSJed Brown           *slot = bcol;
756cd093f1dSJed Brown           seen[bcol] = 1;
757cd093f1dSJed Brown           packlen++;
758cd093f1dSJed Brown         }
759cd093f1dSJed Brown       }
760cd093f1dSJed Brown     }
761cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
762cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
763cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
764cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
765cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
766cd093f1dSJed Brown   }
767cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
768cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
769cd093f1dSJed Brown 
770cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
771cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
772cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
773cd093f1dSJed Brown 
774cd093f1dSJed Brown   /* put together the new symbolic matrix */
775cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
776cd093f1dSJed Brown 
777cd093f1dSJed Brown   (*C)->rmap->bs = A->rmap->bs;
778cd093f1dSJed Brown   (*C)->cmap->bs = B->cmap->bs;
779cd093f1dSJed Brown 
780cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
781cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
782cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
783cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
784cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
785cd093f1dSJed Brown   c->nonew   = 0;
786cd093f1dSJed Brown 
787cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
788cd093f1dSJed Brown 
789cd093f1dSJed Brown   /* set MatInfo */
790cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
791cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
792cd093f1dSJed Brown   c->maxnz                     = ci[am];
793cd093f1dSJed Brown   c->nz                        = ci[am];
794cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
795cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
796cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
797cd093f1dSJed Brown 
798cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
799cd093f1dSJed Brown   if (ci[am]) {
800cd093f1dSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
801cd093f1dSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
802cd093f1dSJed Brown   } else {
803cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
804cd093f1dSJed Brown   }
805cd093f1dSJed Brown #endif
806cd093f1dSJed Brown   PetscFunctionReturn(0);
807cd093f1dSJed Brown }
808cd093f1dSJed Brown 
809d2853540SHong Zhang /* This routine is not used. Should be removed! */
810bc011b1eSHong Zhang #undef __FUNCT__
8116fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
8126fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8135df89d91SHong Zhang {
814bc011b1eSHong Zhang   PetscErrorCode ierr;
815bc011b1eSHong Zhang 
816bc011b1eSHong Zhang   PetscFunctionBegin;
817bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
8183ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
8196fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8203ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
821bc011b1eSHong Zhang   }
8223ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
8236fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8243ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
825bc011b1eSHong Zhang   PetscFunctionReturn(0);
826bc011b1eSHong Zhang }
827bc011b1eSHong Zhang 
828bc011b1eSHong Zhang #undef __FUNCT__
829e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
830e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
8312128a86cSHong Zhang {
8322128a86cSHong Zhang   PetscErrorCode      ierr;
833e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
8342128a86cSHong Zhang 
8352128a86cSHong Zhang   PetscFunctionBegin;
836b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
8372128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
8382128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
8392128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
8402128a86cSHong Zhang   PetscFunctionReturn(0);
8412128a86cSHong Zhang }
8422128a86cSHong Zhang 
8432128a86cSHong Zhang #undef __FUNCT__
8442128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
8452128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8462128a86cSHong Zhang {
8472128a86cSHong Zhang   PetscErrorCode      ierr;
8482128a86cSHong Zhang   PetscContainer      container;
8490298fd71SBarry Smith   Mat_MatMatTransMult *multtrans=NULL;
8502128a86cSHong Zhang 
8512128a86cSHong Zhang   PetscFunctionBegin;
852e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr);
8532128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
8542128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr);
8552205254eSKarl Rupp 
8562128a86cSHong Zhang   A->ops->destroy = multtrans->destroy;
8572128a86cSHong Zhang   if (A->ops->destroy) {
8582128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
8592128a86cSHong Zhang   }
860e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
8612128a86cSHong Zhang   PetscFunctionReturn(0);
8622128a86cSHong Zhang }
8632128a86cSHong Zhang 
8642128a86cSHong Zhang #undef __FUNCT__
8656fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
8666fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
867bc011b1eSHong Zhang {
8685df89d91SHong Zhang   PetscErrorCode      ierr;
86981d82fe4SHong Zhang   Mat                 Bt;
87081d82fe4SHong Zhang   PetscInt            *bti,*btj;
871e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
8722128a86cSHong Zhang   PetscContainer      container;
873d2853540SHong Zhang 
87481d82fe4SHong Zhang   PetscFunctionBegin;
87581d82fe4SHong Zhang   /* create symbolic Bt */
87681d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8770298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
8782205254eSKarl Rupp 
879c3e5699bSSatish Balay   Bt->rmap->bs = A->cmap->bs;
880c3e5699bSSatish Balay   Bt->cmap->bs = B->cmap->bs;
88181d82fe4SHong Zhang 
88281d82fe4SHong Zhang   /* get symbolic C=A*Bt */
88381d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
88481d82fe4SHong Zhang 
8852128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
886e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
8872128a86cSHong Zhang 
8882128a86cSHong Zhang   /* attach the supporting struct to C */
8892128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
8902128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
891e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
892e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
8932128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
8942128a86cSHong Zhang 
8952128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
8962128a86cSHong Zhang   multtrans->destroy     = (*C)->ops->destroy;
8972128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8982128a86cSHong Zhang 
8990298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-matmattransmult_color",&multtrans->usecoloring,NULL);CHKERRQ(ierr);
9002128a86cSHong Zhang   if (multtrans->usecoloring) {
901b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
902b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
9032128a86cSHong Zhang     ISColoring           iscoloring;
9042128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
905*4d478ae7SHong Zhang     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
906*4d478ae7SHong Zhang     /* inode causes memory problem, don't know why */
907*4d478ae7SHong Zhang     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
908e8354b3eSHong Zhang 
9094614e006SHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
910b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
9112205254eSKarl Rupp 
9122128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
9132205254eSKarl Rupp 
9142128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
9152128a86cSHong Zhang 
9162128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
9172128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
9182128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9192128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
9200298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9212205254eSKarl Rupp 
9222128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
9232128a86cSHong Zhang     multtrans->Bt_den   = Bt_dense;
9242128a86cSHong Zhang 
9252128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9262128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9272128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9280298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9292205254eSKarl Rupp 
9302128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
9312128a86cSHong Zhang     multtrans->ABt_den  = C_dense;
932f94ccd6cSHong Zhang 
933f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9341ce71dffSSatish Balay     {
935f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
93622d28d08SBarry Smith       ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr);
9371ce71dffSSatish Balay     }
938f94ccd6cSHong Zhang #endif
9392128a86cSHong Zhang   }
940e99be685SHong Zhang   /* clean up */
941e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
942e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9435df89d91SHong Zhang   PetscFunctionReturn(0);
9445df89d91SHong Zhang }
9455df89d91SHong Zhang 
9466973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
9475df89d91SHong Zhang #undef __FUNCT__
9486fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9496fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9505df89d91SHong Zhang {
9515df89d91SHong Zhang   PetscErrorCode      ierr;
9525df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
953e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
95481d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9555df89d91SHong Zhang   PetscLogDouble      flops=0.0;
956aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
957e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
9582128a86cSHong Zhang   PetscContainer      container;
9596973769bSHong Zhang #if defined(USE_ARRAY)
9606973769bSHong Zhang   MatScalar *spdot;
9616973769bSHong Zhang #endif
9625df89d91SHong Zhang 
9635df89d91SHong Zhang   PetscFunctionBegin;
96458ed3ceaSHong Zhang   /* clear old values in C */
96558ed3ceaSHong Zhang   if (!c->a) {
96658ed3ceaSHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
96758ed3ceaSHong Zhang     c->a      = ca;
96858ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
96958ed3ceaSHong Zhang   } else {
97058ed3ceaSHong Zhang     ca =  c->a;
97158ed3ceaSHong Zhang   }
97258ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
97358ed3ceaSHong Zhang 
974e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr);
9752128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
9762128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr);
9772128a86cSHong Zhang   if (multtrans->usecoloring) {
978b9af6bddSHong Zhang     MatTransposeColoring matcoloring = multtrans->matcoloring;
979c8db22f6SHong Zhang     Mat                  Bt_dense;
980c8db22f6SHong Zhang     PetscInt             m,n;
981a0b95ffbSSatish Balay     Mat                  C_dense = multtrans->ABt_den;
982b3b4fc5aSHong Zhang     PetscLogDouble       t0,t1,t2,t3,Bt_den,C_den,C_sp;
983c8db22f6SHong Zhang 
984b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
98526be7272SHong Zhang     ierr = PetscTime(&t0);CHKERRQ(ierr);
986b3b4fc5aSHong Zhang     Bt_dense = multtrans->Bt_den;
987b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
98826be7272SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
989b3b4fc5aSHong Zhang     Bt_den = t1 - t0;
990c8db22f6SHong Zhang 
991c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9922128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
99326be7272SHong Zhang     ierr = PetscTime(&t2);CHKERRQ(ierr);
994b3b4fc5aSHong Zhang     C_den = t2 - t1;
995c8db22f6SHong Zhang 
9962128a86cSHong Zhang     /* Recover C from C_dense */
997b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
99826be7272SHong Zhang     ierr = PetscTime(&t3);CHKERRQ(ierr);
999b3b4fc5aSHong Zhang     C_sp = t3 - t2;
1000b3b4fc5aSHong Zhang #if defined(PETSC_USE_INFO)
1001b3b4fc5aSHong Zhang     ierr = PetscInfo4(C,"Coloring A*B^T = Bt_den %g + C_den %g + C_sp %g = %g;\n",Bt_den,C_den,C_sp,Bt_den+C_den+C_sp);CHKERRQ(ierr);
1002b3b4fc5aSHong Zhang     ierr     = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
1003b3b4fc5aSHong Zhang     ierr = PetscInfo4(C,"Bt_den: %d x %d, B: %d x %d\n",m,n,m,C->cmap->n);CHKERRQ(ierr);
1004b3b4fc5aSHong Zhang #endif
1005c8db22f6SHong Zhang     PetscFunctionReturn(0);
1006c8db22f6SHong Zhang   }
1007c8db22f6SHong Zhang 
10086973769bSHong Zhang #if defined(USE_ARRAY)
10096973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
1010e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
1011e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
10126973769bSHong Zhang #endif
10136973769bSHong Zhang 
10141fa1209cSHong Zhang   for (i=0; i<cm; i++) {
101581d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
101681d82fe4SHong Zhang     acol = aj + ai[i];
10176973769bSHong Zhang     aval = aa + ai[i];
10181fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
10191fa1209cSHong Zhang     ccol = cj + ci[i];
10206973769bSHong Zhang     cval = ca + ci[i];
10211fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
102281d82fe4SHong Zhang       brow = ccol[j];
102381d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
102481d82fe4SHong Zhang       bcol = bj + bi[brow];
10256973769bSHong Zhang       bval = ba + bi[brow];
10266973769bSHong Zhang 
102781d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
10286973769bSHong Zhang #if defined(USE_ARRAY)
10296973769bSHong Zhang       /* put ba to spdot array */
10306973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
10316973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
10326973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++) {
10336973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
10346973769bSHong Zhang       }
10356973769bSHong Zhang       /* zero spdot array */
10366973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
10376973769bSHong Zhang #else
103881d82fe4SHong Zhang       nexta = 0; nextb = 0;
103981d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
10407b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
104181d82fe4SHong Zhang         if (nexta == anzi) break;
10427b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
104381d82fe4SHong Zhang         if (nextb == bnzj) break;
104481d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
10456973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
104681d82fe4SHong Zhang           nexta++; nextb++;
104781d82fe4SHong Zhang           flops += 2;
10481fa1209cSHong Zhang         }
10491fa1209cSHong Zhang       }
10506973769bSHong Zhang #endif
105181d82fe4SHong Zhang     }
105281d82fe4SHong Zhang   }
105381d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105481d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105581d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10566973769bSHong Zhang #if defined(USE_ARRAY)
10576973769bSHong Zhang   ierr = PetscFree(spdot);
10586973769bSHong Zhang #endif
10595df89d91SHong Zhang   PetscFunctionReturn(0);
10605df89d91SHong Zhang }
10615df89d91SHong Zhang 
10625df89d91SHong Zhang #undef __FUNCT__
106375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
10640adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10650adebc6cSBarry Smith {
10665df89d91SHong Zhang   PetscErrorCode ierr;
10675df89d91SHong Zhang 
10685df89d91SHong Zhang   PetscFunctionBegin;
10695df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
107075648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
10715df89d91SHong Zhang   }
107275648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
10735df89d91SHong Zhang   PetscFunctionReturn(0);
10745df89d91SHong Zhang }
10755df89d91SHong Zhang 
10765df89d91SHong Zhang #undef __FUNCT__
107775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
107875648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10795df89d91SHong Zhang {
1080bc011b1eSHong Zhang   PetscErrorCode ierr;
1081bc011b1eSHong Zhang   Mat            At;
108238baddfdSBarry Smith   PetscInt       *ati,*atj;
1083bc011b1eSHong Zhang 
1084bc011b1eSHong Zhang   PetscFunctionBegin;
1085bc011b1eSHong Zhang   /* create symbolic At */
1086bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10870298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
10882205254eSKarl Rupp 
1089a2f3521dSMark F. Adams   At->rmap->bs = A->cmap->bs;
1090a2f3521dSMark F. Adams   At->cmap->bs = B->cmap->bs;
1091bc011b1eSHong Zhang 
1092bc011b1eSHong Zhang   /* get symbolic C=At*B */
1093bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1094bc011b1eSHong Zhang 
1095bc011b1eSHong Zhang   /* clean up */
10966bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1097bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1098bc011b1eSHong Zhang   PetscFunctionReturn(0);
1099bc011b1eSHong Zhang }
1100bc011b1eSHong Zhang 
1101bc011b1eSHong Zhang #undef __FUNCT__
110275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
110375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1104bc011b1eSHong Zhang {
1105bc011b1eSHong Zhang   PetscErrorCode ierr;
11060fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1107d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1108d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1109d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1110aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1111bc011b1eSHong Zhang 
1112bc011b1eSHong Zhang   PetscFunctionBegin;
1113aa1aec99SHong Zhang   if (!c->a) {
1114aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
11152205254eSKarl Rupp 
1116aa1aec99SHong Zhang     c->a      = ca;
1117aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1118aa1aec99SHong Zhang   } else {
1119aa1aec99SHong Zhang     ca = c->a;
1120aa1aec99SHong Zhang   }
1121bc011b1eSHong Zhang   /* clear old values in C */
1122bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1123bc011b1eSHong Zhang 
1124bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1125bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1126bc011b1eSHong Zhang     bj   = b->j + bi[i];
1127bc011b1eSHong Zhang     ba   = b->a + bi[i];
1128bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1129bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1130bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1131bc011b1eSHong Zhang       nextb = 0;
11320fbc74f4SHong Zhang       crow  = *aj++;
1133bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1134bc011b1eSHong Zhang       caj   = ca + ci[crow];
1135bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1136bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
11370fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
11380fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1139bc011b1eSHong Zhang           nextb++;
1140bc011b1eSHong Zhang         }
1141bc011b1eSHong Zhang       }
1142bc011b1eSHong Zhang       flops += 2*bnzi;
11430fbc74f4SHong Zhang       aa++;
1144bc011b1eSHong Zhang     }
1145bc011b1eSHong Zhang   }
1146bc011b1eSHong Zhang 
1147bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1148bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1149bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1150bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1151bc011b1eSHong Zhang   PetscFunctionReturn(0);
1152bc011b1eSHong Zhang }
11539123193aSHong Zhang 
11549123193aSHong Zhang #undef __FUNCT__
11559123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
11569123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11579123193aSHong Zhang {
11589123193aSHong Zhang   PetscErrorCode ierr;
11599123193aSHong Zhang 
11609123193aSHong Zhang   PetscFunctionBegin;
11619123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
11623ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11639123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11643ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11659123193aSHong Zhang   }
11663ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11679123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11684614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11699123193aSHong Zhang   PetscFunctionReturn(0);
11709123193aSHong Zhang }
1171edd81eecSMatthew Knepley 
11729123193aSHong Zhang #undef __FUNCT__
11739123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11749123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11759123193aSHong Zhang {
11769123193aSHong Zhang   PetscErrorCode ierr;
11779123193aSHong Zhang 
11789123193aSHong Zhang   PetscFunctionBegin;
11795a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11802205254eSKarl Rupp 
1181d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11829123193aSHong Zhang   PetscFunctionReturn(0);
11839123193aSHong Zhang }
11849123193aSHong Zhang 
11859123193aSHong Zhang #undef __FUNCT__
11869123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11879123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11889123193aSHong Zhang {
1189f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11909123193aSHong Zhang   PetscErrorCode ierr;
1191dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1192dd6ea824SBarry Smith   MatScalar      *aa;
1193d0f46423SBarry Smith   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
1194f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
11959123193aSHong Zhang 
11969123193aSHong Zhang   PetscFunctionBegin;
1197f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1198e32f2f54SBarry Smith   if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm);
1199e32f2f54SBarry 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);
1200e32f2f54SBarry 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);
12018c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12028c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1203f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1204f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1205f32d5d43SBarry Smith     colam = col*am;
1206f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1207f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1208f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1209f32d5d43SBarry Smith       aj = a->j + a->i[i];
1210f32d5d43SBarry Smith       aa = a->a + a->i[i];
1211f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1212f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
1213f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
1214f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
1215f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
12169123193aSHong Zhang       }
1217f32d5d43SBarry Smith       c[colam + i]       = r1;
1218f32d5d43SBarry Smith       c[colam + am + i]  = r2;
1219f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
1220f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
1221f32d5d43SBarry Smith     }
1222f32d5d43SBarry Smith     b1 += bm4;
1223f32d5d43SBarry Smith     b2 += bm4;
1224f32d5d43SBarry Smith     b3 += bm4;
1225f32d5d43SBarry Smith     b4 += bm4;
1226f32d5d43SBarry Smith   }
1227f32d5d43SBarry Smith   for (; col<cn; col++) {     /* over extra columns of C */
1228f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1229f32d5d43SBarry Smith       r1 = 0.0;
1230f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1231f32d5d43SBarry Smith       aj = a->j + a->i[i];
1232f32d5d43SBarry Smith       aa = a->a + a->i[i];
1233f32d5d43SBarry Smith 
1234f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1235f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
1236f32d5d43SBarry Smith       }
1237f32d5d43SBarry Smith       c[col*am + i] = r1;
1238f32d5d43SBarry Smith     }
1239f32d5d43SBarry Smith     b1 += bm;
1240f32d5d43SBarry Smith   }
1241dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
12428c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
12438c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1244f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1245f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1246f32d5d43SBarry Smith   PetscFunctionReturn(0);
1247f32d5d43SBarry Smith }
1248f32d5d43SBarry Smith 
1249f32d5d43SBarry Smith /*
12504324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1251f32d5d43SBarry Smith */
1252f32d5d43SBarry Smith #undef __FUNCT__
1253f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1254f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1255f32d5d43SBarry Smith {
1256f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1257f32d5d43SBarry Smith   PetscErrorCode ierr;
1258dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1259dd6ea824SBarry Smith   MatScalar      *aa;
1260d0f46423SBarry Smith   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
12614324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1262f32d5d43SBarry Smith 
1263f32d5d43SBarry Smith   PetscFunctionBegin;
1264f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
12658c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12668c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1267f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12684324174eSBarry Smith 
12694324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12704324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12714324174eSBarry Smith       colam = col*am;
12724324174eSBarry Smith       arm   = a->compressedrow.nrows;
12734324174eSBarry Smith       ii    = a->compressedrow.i;
12744324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12754324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12764324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12774324174eSBarry Smith         n  = ii[i+1] - ii[i];
12784324174eSBarry Smith         aj = a->j + ii[i];
12794324174eSBarry Smith         aa = a->a + ii[i];
12804324174eSBarry Smith         for (j=0; j<n; j++) {
12814324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12824324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12834324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12844324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12854324174eSBarry Smith         }
12864324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12874324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12884324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12894324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12904324174eSBarry Smith       }
12914324174eSBarry Smith       b1 += bm4;
12924324174eSBarry Smith       b2 += bm4;
12934324174eSBarry Smith       b3 += bm4;
12944324174eSBarry Smith       b4 += bm4;
12954324174eSBarry Smith     }
12964324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12974324174eSBarry Smith       colam = col*am;
12984324174eSBarry Smith       arm   = a->compressedrow.nrows;
12994324174eSBarry Smith       ii    = a->compressedrow.i;
13004324174eSBarry Smith       ridx  = a->compressedrow.rindex;
13014324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
13024324174eSBarry Smith         r1 = 0.0;
13034324174eSBarry Smith         n  = ii[i+1] - ii[i];
13044324174eSBarry Smith         aj = a->j + ii[i];
13054324174eSBarry Smith         aa = a->a + ii[i];
13064324174eSBarry Smith 
13074324174eSBarry Smith         for (j=0; j<n; j++) {
13084324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
13094324174eSBarry Smith         }
1310a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
13114324174eSBarry Smith       }
13124324174eSBarry Smith       b1 += bm;
13134324174eSBarry Smith     }
13144324174eSBarry Smith   } else {
1315f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1316f32d5d43SBarry Smith       colam = col*am;
1317f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1318f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1319f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1320f32d5d43SBarry Smith         aj = a->j + a->i[i];
1321f32d5d43SBarry Smith         aa = a->a + a->i[i];
1322f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1323f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1324f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1325f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1326f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1327f32d5d43SBarry Smith         }
1328f32d5d43SBarry Smith         c[colam + i]       += r1;
1329f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1330f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1331f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1332f32d5d43SBarry Smith       }
1333f32d5d43SBarry Smith       b1 += bm4;
1334f32d5d43SBarry Smith       b2 += bm4;
1335f32d5d43SBarry Smith       b3 += bm4;
1336f32d5d43SBarry Smith       b4 += bm4;
1337f32d5d43SBarry Smith     }
1338f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1339a2ea699eSBarry Smith       colam = col*am;
1340f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1341f32d5d43SBarry Smith         r1 = 0.0;
1342f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1343f32d5d43SBarry Smith         aj = a->j + a->i[i];
1344f32d5d43SBarry Smith         aa = a->a + a->i[i];
1345f32d5d43SBarry Smith 
1346f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1347f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1348f32d5d43SBarry Smith         }
1349a2ea699eSBarry Smith         c[colam + i] += r1;
1350f32d5d43SBarry Smith       }
1351f32d5d43SBarry Smith       b1 += bm;
1352f32d5d43SBarry Smith     }
13534324174eSBarry Smith   }
1354dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13558c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
13568c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13579123193aSHong Zhang   PetscFunctionReturn(0);
13589123193aSHong Zhang }
1359b1683b59SHong Zhang 
1360b1683b59SHong Zhang #undef __FUNCT__
1361b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1362b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1363c8db22f6SHong Zhang {
1364c8db22f6SHong Zhang   PetscErrorCode ierr;
13652f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13662f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13672f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13682f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13692f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13702f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1371c8db22f6SHong Zhang 
1372c8db22f6SHong Zhang   PetscFunctionBegin;
13732f5f1f90SHong Zhang   btval_den=btdense->v;
13742f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13752f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13762f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13772f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1378d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13792f5f1f90SHong Zhang       btcol = bj + bi[col];
13802f5f1f90SHong Zhang       btval = ba + bi[col];
13812f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1382d2853540SHong Zhang       for (j=0; j<anz; j++) {
13832f5f1f90SHong Zhang         brow            = btcol[j];
13842f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1385c8db22f6SHong Zhang       }
1386c8db22f6SHong Zhang     }
13872f5f1f90SHong Zhang     btval_den += m;
1388c8db22f6SHong Zhang   }
1389c8db22f6SHong Zhang   PetscFunctionReturn(0);
1390c8db22f6SHong Zhang }
1391c8db22f6SHong Zhang 
1392c8db22f6SHong Zhang #undef __FUNCT__
1393b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1394b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13958972f759SHong Zhang {
13968972f759SHong Zhang   PetscErrorCode ierr;
1397b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
13982f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1399b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
14002f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
14018972f759SHong Zhang 
14028972f759SHong Zhang   PetscFunctionBegin;
14030298fd71SBarry Smith   ierr   = MatGetLocalSize(Csp,&m,NULL);CHKERRQ(ierr);
1404a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1405b2d2b04fSHong Zhang   cp_den = ca_den;
14062f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
14072f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
14082f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
14092f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
14102f5f1f90SHong Zhang     for (l=0; l<nrows; l++) {
14112f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1412b2d2b04fSHong Zhang     }
1413b2d2b04fSHong Zhang     cp_den += m;
1414b2d2b04fSHong Zhang   }
1415a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
14168972f759SHong Zhang   PetscFunctionReturn(0);
14178972f759SHong Zhang }
14188972f759SHong Zhang 
1419e99be685SHong Zhang /*
1420e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1421e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1422e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1423e99be685SHong Zhang  */
1424e99be685SHong Zhang #undef __FUNCT__
1425e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
14261a83f524SJed Brown PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1427e99be685SHong Zhang {
1428e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1429e99be685SHong Zhang   PetscErrorCode ierr;
1430e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1431e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1432e99be685SHong Zhang   PetscInt       *cspidx;
1433e99be685SHong Zhang 
1434e99be685SHong Zhang   PetscFunctionBegin;
1435e99be685SHong Zhang   *nn = n;
1436e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1437e99be685SHong Zhang   if (symmetric) {
1438ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
14391a83f524SJed Brown     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr);
1440e99be685SHong Zhang   } else {
1441e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1442e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1443e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1444e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1445e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1446e99be685SHong Zhang     jj   = a->j;
1447e99be685SHong Zhang     for (i=0; i<nz; i++) {
1448e99be685SHong Zhang       collengths[jj[i]]++;
1449e99be685SHong Zhang     }
1450e99be685SHong Zhang     cia[0] = oshift;
1451e99be685SHong Zhang     for (i=0; i<n; i++) {
1452e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1453e99be685SHong Zhang     }
1454e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1455e99be685SHong Zhang     jj   = a->j;
1456e99be685SHong Zhang     for (row=0; row<m; row++) {
1457e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1458e99be685SHong Zhang       for (i=0; i<mr; i++) {
1459e99be685SHong Zhang         col = *jj++;
14602205254eSKarl Rupp 
1461e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1462e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift]  = row + oshift;
1463e99be685SHong Zhang       }
1464e99be685SHong Zhang     }
1465e99be685SHong Zhang     ierr   = PetscFree(collengths);CHKERRQ(ierr);
1466e99be685SHong Zhang     *ia    = cia; *ja = cja;
1467e99be685SHong Zhang     *spidx = cspidx;
1468e99be685SHong Zhang   }
1469e99be685SHong Zhang   PetscFunctionReturn(0);
1470e99be685SHong Zhang }
1471e99be685SHong Zhang 
1472e99be685SHong Zhang #undef __FUNCT__
1473e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
14741a83f524SJed Brown PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1475e99be685SHong Zhang {
1476e99be685SHong Zhang   PetscErrorCode ierr;
1477e99be685SHong Zhang 
1478e99be685SHong Zhang   PetscFunctionBegin;
1479e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1480e99be685SHong Zhang 
1481e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1482e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1483e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1484e99be685SHong Zhang   PetscFunctionReturn(0);
1485e99be685SHong Zhang }
1486e99be685SHong Zhang 
14878972f759SHong Zhang #undef __FUNCT__
1488b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1489b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1490b1683b59SHong Zhang {
1491b1683b59SHong Zhang   PetscErrorCode ierr;
14921a83f524SJed Brown   PetscInt       i,n,nrows,N,j,k,m,ncols,col,cm;
14931a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1494b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1495b1683b59SHong Zhang   IS             *isa;
1496b1683b59SHong Zhang   PetscBool      flg1,flg2;
1497b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1498e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1499d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1500b1683b59SHong Zhang 
1501b1683b59SHong Zhang   PetscFunctionBegin;
1502b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1503e99be685SHong Zhang 
1504b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1505251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1506251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1507b1683b59SHong Zhang   if (flg1 || flg2) {
1508b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1509b1683b59SHong Zhang   }
1510b1683b59SHong Zhang 
1511b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1512b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1513b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1514b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1515b1683b59SHong Zhang   c->rstart = 0;
1516b1683b59SHong Zhang 
1517b1683b59SHong Zhang   c->ncolors = nis;
1518b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1519b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1520d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1521d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
15222205254eSKarl Rupp 
1523d2853540SHong Zhang   colorforrow[0]    = 0;
1524d2853540SHong Zhang   rows_i            = rows;
1525d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1526b1683b59SHong Zhang 
1527d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1528d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
15292205254eSKarl Rupp 
1530d2853540SHong Zhang   colorforcol[0] = 0;
1531d2853540SHong Zhang   columns_i      = columns;
1532d2853540SHong Zhang 
15330298fd71SBarry Smith   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); /* column-wise storage of mat */
1534b1683b59SHong Zhang 
1535b28d80bdSHong Zhang   cm   = c->m;
1536b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1537e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1538b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1539b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1540b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
15412205254eSKarl Rupp 
1542b1683b59SHong Zhang     c->ncolumns[i] = n;
1543b1683b59SHong Zhang     if (n) {
1544d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1545b1683b59SHong Zhang     }
1546d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1547d2853540SHong Zhang     columns_i       += n;
1548b1683b59SHong Zhang 
1549b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1550e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1551e99be685SHong Zhang 
1552b1683b59SHong Zhang     /* loop over columns*/
1553b1683b59SHong Zhang     for (j=0; j<n; j++) {
1554b1683b59SHong Zhang       col     = is[j];
1555d2853540SHong Zhang       row_idx = cj + ci[col];
1556b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1557b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1558b1683b59SHong Zhang       for (k=0; k<m; k++) {
1559e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1560d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1561b1683b59SHong Zhang       }
1562b1683b59SHong Zhang     }
1563b1683b59SHong Zhang     /* count the number of hits */
1564b1683b59SHong Zhang     nrows = 0;
1565e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1566b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1567b1683b59SHong Zhang     }
1568b1683b59SHong Zhang     c->nrows[i]      = nrows;
1569d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1570d2853540SHong Zhang 
1571b1683b59SHong Zhang     nrows = 0;
1572e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1573b1683b59SHong Zhang       if (rowhit[j]) {
1574d2853540SHong Zhang         rows_i[nrows]            = j;
1575e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1576b1683b59SHong Zhang         nrows++;
1577b1683b59SHong Zhang       }
1578b1683b59SHong Zhang     }
1579b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1580d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1581b1683b59SHong Zhang   }
15820298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1583b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1584b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1585d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1586d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1587d2853540SHong Zhang #endif
1588b28d80bdSHong Zhang 
1589d2853540SHong Zhang   c->colorforrow     = colorforrow;
1590d2853540SHong Zhang   c->rows            = rows;
1591d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1592d2853540SHong Zhang   c->colorforcol     = colorforcol;
1593d2853540SHong Zhang   c->columns         = columns;
15942205254eSKarl Rupp 
1595f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1596b1683b59SHong Zhang   PetscFunctionReturn(0);
1597b1683b59SHong Zhang }
1598