xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision b3b4fc5a114eeab53373cded8564dfacad3aa869)
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>
127bab7c10SHong Zhang 
1358cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
14cd093f1dSJed Brown 
156284ec50SHong Zhang #undef __FUNCT__
165c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1738baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1838baddfdSBarry Smith {
19dfbe8321SBarry Smith   PetscErrorCode ierr;
2058cf0668SJed Brown   PetscBool      scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE,btheap = PETSC_FALSE,llcondensed = PETSC_FALSE;
215c66b693SKris Buschelman 
225c66b693SKris Buschelman   PetscFunctionBegin;
2326be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
24152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
250298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,NULL);CHKERRQ(ierr);
260298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,NULL);CHKERRQ(ierr);
270298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,NULL);CHKERRQ(ierr);
280298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_btheap","Use btheap implementation of symbolic factorization C=A*B","",btheap,&btheap,NULL);CHKERRQ(ierr);
2958cf0668SJed Brown     ierr = PetscOptionsBool("-matmatmult_llcondensed","Use LLCondensed to for symbolic C=A*B","",llcondensed,&llcondensed,NULL);CHKERRQ(ierr);
30d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
313ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
323c50cad2SHong Zhang     if (scalable_fast) {
333c50cad2SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
34aacf854cSBarry Smith     } else if (scalable) {
35aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
360ced3a2bSJed Brown     } else if (heap) {
370ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
388a07c6f1SJed Brown     } else if (btheap) {
398a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
4058cf0668SJed Brown     } else if (llcondensed) {
4158cf0668SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
4225023028SHong Zhang     } else {
4326be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
4425023028SHong Zhang     }
453ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
4626be0446SHong Zhang   }
475c913ed7SHong Zhang 
483ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
4901e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
503ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
515c66b693SKris Buschelman   PetscFunctionReturn(0);
525c66b693SKris Buschelman }
531c24bd37SHong Zhang 
54e9e4536cSHong Zhang #undef __FUNCT__
5558cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed"
5658cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
57b561aa9dSHong Zhang {
58b561aa9dSHong Zhang   PetscErrorCode     ierr;
59b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
60971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
61c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
62b561aa9dSHong Zhang   PetscReal          afill;
63fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
64fb908031SHong Zhang   PetscBT            lnkbt;
650298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
66b561aa9dSHong Zhang 
67b561aa9dSHong Zhang   PetscFunctionBegin;
68bd958071SHong Zhang   /* Get ci and cj */
69bd958071SHong Zhang   /*---------------*/
70fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
71fb908031SHong Zhang   /* free space for accumulating nonzero column info */
72fb908031SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
73fb908031SHong Zhang   ci[0] = 0;
74fb908031SHong Zhang 
75fb908031SHong Zhang   /* create and initialize a linked list */
76fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
77fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
78fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
79fb908031SHong Zhang 
80fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
81fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
822205254eSKarl Rupp 
83fb908031SHong Zhang   current_space = free_space;
84fb908031SHong Zhang 
85fb908031SHong Zhang   /* Determine ci and cj */
86fb908031SHong Zhang   for (i=0; i<am; i++) {
87fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
88fb908031SHong Zhang     aj   = a->j + ai[i];
89fb908031SHong Zhang     for (j=0; j<anzi; j++) {
90fb908031SHong Zhang       brow = aj[j];
91fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
92fb908031SHong Zhang       bj   = b->j + bi[brow];
93fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
94fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
95fb908031SHong Zhang     }
96fb908031SHong Zhang     cnzi = lnk[0];
97fb908031SHong Zhang 
98fb908031SHong Zhang     /* If free space is not available, make more free space */
99fb908031SHong Zhang     /* Double the amount of total space in the list */
100fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
101fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
102fb908031SHong Zhang       ndouble++;
103fb908031SHong Zhang     }
104fb908031SHong Zhang 
105fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
106fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1072205254eSKarl Rupp 
108fb908031SHong Zhang     current_space->array           += cnzi;
109fb908031SHong Zhang     current_space->local_used      += cnzi;
110fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1112205254eSKarl Rupp 
112fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
113fb908031SHong Zhang   }
114fb908031SHong Zhang 
115fb908031SHong Zhang   /* Column indices are in the list of free space */
116fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
117fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
118fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
119fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
120fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
121b561aa9dSHong Zhang 
12226be0446SHong Zhang   /* put together the new symbolic matrix */
123ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1242205254eSKarl Rupp 
125a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
126a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
12758c24d83SHong Zhang 
12858c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
12958c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
13058c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
131aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
132e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
13358c24d83SHong Zhang   c->nonew                  = 0;
134002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1350b7e3e3dSHong Zhang 
1367212da7cSHong Zhang   /* set MatInfo */
1377212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
138f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1397212da7cSHong Zhang   c->maxnz                     = ci[am];
1407212da7cSHong Zhang   c->nz                        = ci[am];
141fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1427212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1437212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1447212da7cSHong Zhang 
1457212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1467212da7cSHong Zhang   if (ci[am]) {
147fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
148f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
149f2b054eeSHong Zhang   } else {
150f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
151be0fcf8dSHong Zhang   }
152f2b054eeSHong Zhang #endif
15358c24d83SHong Zhang   PetscFunctionReturn(0);
15458c24d83SHong Zhang }
155d50806bdSBarry Smith 
15626be0446SHong Zhang #undef __FUNCT__
15726be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
158dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
159d50806bdSBarry Smith {
160dfbe8321SBarry Smith   PetscErrorCode ierr;
161d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
162d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
163d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
164d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
16538baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
166c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
167519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
168aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
169aa1aec99SHong Zhang   PetscScalar    *ab_dense;
170d50806bdSBarry Smith 
171d50806bdSBarry Smith   PetscFunctionBegin;
172aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
173aa1aec99SHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
174aa1aec99SHong Zhang     c->a      = ca;
175aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
176aa1aec99SHong Zhang 
177aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
178aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
1792205254eSKarl Rupp 
180aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
181aa1aec99SHong Zhang   } else {
182aa1aec99SHong Zhang     ca       = c->a;
183aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
184aa1aec99SHong Zhang   }
185aa1aec99SHong Zhang 
186c124e916SHong Zhang   /* clean old values in C */
187c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
188d50806bdSBarry Smith   /* Traverse A row-wise. */
189d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
190d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
191d50806bdSBarry Smith   for (i=0; i<am; i++) {
192d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
193d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
194519eb980SHong Zhang       brow = aj[j];
195d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
196d50806bdSBarry Smith       bjj  = bj + bi[brow];
197d50806bdSBarry Smith       baj  = ba + bi[brow];
198519eb980SHong Zhang       /* perform dense axpy */
19936ec6d2dSHong Zhang       valtmp = aa[j];
200519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
20136ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
202519eb980SHong Zhang       }
203519eb980SHong Zhang       flops += 2*bnzi;
204519eb980SHong Zhang     }
205c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
206c58ca2e3SHong Zhang 
207c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
208519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
209519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
210519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
211519eb980SHong Zhang     }
212519eb980SHong Zhang     flops += cnzi;
213519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
214519eb980SHong Zhang   }
215c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
216c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
217c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
218c58ca2e3SHong Zhang   PetscFunctionReturn(0);
219c58ca2e3SHong Zhang }
220c58ca2e3SHong Zhang 
221c58ca2e3SHong Zhang #undef __FUNCT__
22225023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
22325023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
224c58ca2e3SHong Zhang {
225c58ca2e3SHong Zhang   PetscErrorCode ierr;
226c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
227c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
228c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
229c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
230c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
231c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
232c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
23336ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
234c58ca2e3SHong Zhang   PetscInt       nextb;
235c58ca2e3SHong Zhang 
236c58ca2e3SHong Zhang   PetscFunctionBegin;
237c58ca2e3SHong Zhang   /* clean old values in C */
238c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
239c58ca2e3SHong Zhang   /* Traverse A row-wise. */
240c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
241c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
242519eb980SHong Zhang   for (i=0; i<am; i++) {
243519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
244519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
245519eb980SHong Zhang     for (j=0; j<anzi; j++) {
246519eb980SHong Zhang       brow = aj[j];
247519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
248519eb980SHong Zhang       bjj  = bj + bi[brow];
249519eb980SHong Zhang       baj  = ba + bi[brow];
250519eb980SHong Zhang       /* perform sparse axpy */
25136ec6d2dSHong Zhang       valtmp = aa[j];
252c124e916SHong Zhang       nextb  = 0;
253c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
254c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
25536ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
256c124e916SHong Zhang         }
257d50806bdSBarry Smith       }
258d50806bdSBarry Smith       flops += 2*bnzi;
259d50806bdSBarry Smith     }
260519eb980SHong Zhang     aj += anzi; aa += anzi;
261519eb980SHong Zhang     cj += cnzi; ca += cnzi;
262d50806bdSBarry Smith   }
263c58ca2e3SHong Zhang 
264716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
265716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
266d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
267d50806bdSBarry Smith   PetscFunctionReturn(0);
268d50806bdSBarry Smith }
269bc011b1eSHong Zhang 
270e9e4536cSHong Zhang #undef __FUNCT__
2713c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2723c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
27325296bd5SBarry Smith {
27425296bd5SBarry Smith   PetscErrorCode     ierr;
27525296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
27625296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2773c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
27825296bd5SBarry Smith   MatScalar          *ca;
27925296bd5SBarry Smith   PetscReal          afill;
2803c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2810298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
28225296bd5SBarry Smith 
28325296bd5SBarry Smith   PetscFunctionBegin;
2843c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2853c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2863c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2873c50cad2SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2883c50cad2SHong Zhang   ci[0] = 0;
2893c50cad2SHong Zhang 
2903c50cad2SHong Zhang   /* create and initialize a linked list */
2913c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2923c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2933c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
2943c50cad2SHong Zhang 
2953c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
2963c50cad2SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
2973c50cad2SHong Zhang   current_space = free_space;
2983c50cad2SHong Zhang 
2993c50cad2SHong Zhang   /* Determine ci and cj */
3003c50cad2SHong Zhang   for (i=0; i<am; i++) {
3013c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3023c50cad2SHong Zhang     aj   = a->j + ai[i];
3033c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3043c50cad2SHong Zhang       brow = aj[j];
3053c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3063c50cad2SHong Zhang       bj   = b->j + bi[brow];
3073c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3083c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3093c50cad2SHong Zhang     }
3103c50cad2SHong Zhang     cnzi = lnk[1];
3113c50cad2SHong Zhang 
3123c50cad2SHong Zhang     /* If free space is not available, make more free space */
3133c50cad2SHong Zhang     /* Double the amount of total space in the list */
3143c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3153c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3163c50cad2SHong Zhang       ndouble++;
3173c50cad2SHong Zhang     }
3183c50cad2SHong Zhang 
3193c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3203c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3212205254eSKarl Rupp 
3223c50cad2SHong Zhang     current_space->array           += cnzi;
3233c50cad2SHong Zhang     current_space->local_used      += cnzi;
3243c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3252205254eSKarl Rupp 
3263c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3273c50cad2SHong Zhang   }
3283c50cad2SHong Zhang 
3293c50cad2SHong Zhang   /* Column indices are in the list of free space */
3303c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3313c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3323c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3333c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3343c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
33525296bd5SBarry Smith 
33625296bd5SBarry Smith   /* Allocate space for ca */
33725296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
33825296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
33925296bd5SBarry Smith 
34025296bd5SBarry Smith   /* put together the new symbolic matrix */
341ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
3422205254eSKarl Rupp 
343a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
344a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
34525296bd5SBarry Smith 
34625296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
34725296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
34825296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
34925296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
35025296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
35125296bd5SBarry Smith   c->nonew   = 0;
3522205254eSKarl Rupp 
35325296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
35425296bd5SBarry Smith 
35525296bd5SBarry Smith   /* set MatInfo */
35625296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
35725296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
35825296bd5SBarry Smith   c->maxnz                     = ci[am];
35925296bd5SBarry Smith   c->nz                        = ci[am];
3603c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
36125296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
36225296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
36325296bd5SBarry Smith 
36425296bd5SBarry Smith #if defined(PETSC_USE_INFO)
36525296bd5SBarry Smith   if (ci[am]) {
3663c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
36725296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
36825296bd5SBarry Smith   } else {
36925296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
37025296bd5SBarry Smith   }
37125296bd5SBarry Smith #endif
37225296bd5SBarry Smith   PetscFunctionReturn(0);
37325296bd5SBarry Smith }
37425296bd5SBarry Smith 
37525296bd5SBarry Smith 
37625296bd5SBarry Smith #undef __FUNCT__
37725023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
37825023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
379e9e4536cSHong Zhang {
380e9e4536cSHong Zhang   PetscErrorCode     ierr;
381e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
382bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
38325c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
384e9e4536cSHong Zhang   MatScalar          *ca;
385e9e4536cSHong Zhang   PetscReal          afill;
386bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
3870298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
388e9e4536cSHong Zhang 
389e9e4536cSHong Zhang   PetscFunctionBegin;
390bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
391bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
392bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
393bd958071SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
394bd958071SHong Zhang   ci[0] = 0;
395bd958071SHong Zhang 
396bd958071SHong Zhang   /* create and initialize a linked list */
397bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
398bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
399bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
400bd958071SHong Zhang 
401bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
402bd958071SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
403bd958071SHong Zhang   current_space = free_space;
404bd958071SHong Zhang 
405bd958071SHong Zhang   /* Determine ci and cj */
406bd958071SHong Zhang   for (i=0; i<am; i++) {
407bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
408bd958071SHong Zhang     aj   = a->j + ai[i];
409bd958071SHong Zhang     for (j=0; j<anzi; j++) {
410bd958071SHong Zhang       brow = aj[j];
411bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
412bd958071SHong Zhang       bj   = b->j + bi[brow];
413bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
414bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
415bd958071SHong Zhang     }
416bd958071SHong Zhang     cnzi = lnk[0];
417bd958071SHong Zhang 
418bd958071SHong Zhang     /* If free space is not available, make more free space */
419bd958071SHong Zhang     /* Double the amount of total space in the list */
420bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
421bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
422bd958071SHong Zhang       ndouble++;
423bd958071SHong Zhang     }
424bd958071SHong Zhang 
425bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
426bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4272205254eSKarl Rupp 
428bd958071SHong Zhang     current_space->array           += cnzi;
429bd958071SHong Zhang     current_space->local_used      += cnzi;
430bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4312205254eSKarl Rupp 
432bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
433bd958071SHong Zhang   }
434bd958071SHong Zhang 
435bd958071SHong Zhang   /* Column indices are in the list of free space */
436bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
437bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
438bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
439bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
440bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
441e9e4536cSHong Zhang 
442e9e4536cSHong Zhang   /* Allocate space for ca */
443bd958071SHong Zhang   /*-----------------------*/
444e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
445e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
446e9e4536cSHong Zhang 
447e9e4536cSHong Zhang   /* put together the new symbolic matrix */
448ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
4492205254eSKarl Rupp 
450a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
451a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
452e9e4536cSHong Zhang 
453e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
454e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
455e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
456e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
457e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
458e9e4536cSHong Zhang   c->nonew   = 0;
4592205254eSKarl Rupp 
46025023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
461e9e4536cSHong Zhang 
462e9e4536cSHong Zhang   /* set MatInfo */
463e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
464e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
465e9e4536cSHong Zhang   c->maxnz                     = ci[am];
466e9e4536cSHong Zhang   c->nz                        = ci[am];
467bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
468e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
469e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
470e9e4536cSHong Zhang 
471e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
472e9e4536cSHong Zhang   if (ci[am]) {
473bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
474e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
475e9e4536cSHong Zhang   } else {
476e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
477e9e4536cSHong Zhang   }
478e9e4536cSHong Zhang #endif
479e9e4536cSHong Zhang   PetscFunctionReturn(0);
480e9e4536cSHong Zhang }
481e9e4536cSHong Zhang 
4820ced3a2bSJed Brown #undef __FUNCT__
4830ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
4840ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
4850ced3a2bSJed Brown {
4860ced3a2bSJed Brown   PetscErrorCode     ierr;
4870ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4880ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
4890ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
4900ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
4910ced3a2bSJed Brown   PetscReal          afill;
4920ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
4930298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
4940ced3a2bSJed Brown   PetscHeap          h;
4950ced3a2bSJed Brown 
4960ced3a2bSJed Brown   PetscFunctionBegin;
497cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
4980ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
4990ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
5000ced3a2bSJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5010ced3a2bSJed Brown   ci[0] = 0;
5020ced3a2bSJed Brown 
5030ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
5040ced3a2bSJed Brown   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5050ced3a2bSJed Brown   current_space = free_space;
5060ced3a2bSJed Brown 
5070ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
5080ced3a2bSJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
5090ced3a2bSJed Brown 
5100ced3a2bSJed Brown   /* Determine ci and cj */
5110ced3a2bSJed Brown   for (i=0; i<am; i++) {
5120ced3a2bSJed 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 */
5130ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5140ced3a2bSJed Brown     ci[i+1] = ci[i];
5150ced3a2bSJed Brown     /* Populate the min heap */
5160ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5170ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5180ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5190ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5200ced3a2bSJed Brown       }
5210ced3a2bSJed Brown     }
5220ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5230ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5240ced3a2bSJed Brown     while (j >= 0) {
5250ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
5260ced3a2bSJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
5270ced3a2bSJed Brown         ndouble++;
5280ced3a2bSJed Brown       }
5290ced3a2bSJed Brown       *(current_space->array++) = col;
5300ced3a2bSJed Brown       current_space->local_used++;
5310ced3a2bSJed Brown       current_space->local_remaining--;
5320ced3a2bSJed Brown       ci[i+1]++;
5330ced3a2bSJed Brown 
5340ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5350ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5360ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5370ced3a2bSJed Brown         PetscInt j2,col2;
5380ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5390ced3a2bSJed Brown         if (col2 != col) break;
5400ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5410ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5420ced3a2bSJed Brown       }
5430ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5440ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5450ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5460ced3a2bSJed Brown     }
5470ced3a2bSJed Brown   }
5480ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5490ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5500ced3a2bSJed Brown 
5510ced3a2bSJed Brown   /* Column indices are in the list of free space */
5520ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5530ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
5540ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5550ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5560ced3a2bSJed Brown 
5570ced3a2bSJed Brown   /* put together the new symbolic matrix */
558ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
5592205254eSKarl Rupp 
5600ced3a2bSJed Brown   (*C)->rmap->bs = A->rmap->bs;
5610ced3a2bSJed Brown   (*C)->cmap->bs = B->cmap->bs;
5620ced3a2bSJed Brown 
5630ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5640ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5650ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5660ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5670ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5680ced3a2bSJed Brown   c->nonew   = 0;
56926fbe8dcSKarl Rupp 
57089d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5710ced3a2bSJed Brown 
5720ced3a2bSJed Brown   /* set MatInfo */
5730ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5740ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5750ced3a2bSJed Brown   c->maxnz                     = ci[am];
5760ced3a2bSJed Brown   c->nz                        = ci[am];
5770ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5780ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5790ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5800ced3a2bSJed Brown 
5810ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5820ced3a2bSJed Brown   if (ci[am]) {
5830ced3a2bSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
5840ced3a2bSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
5850ced3a2bSJed Brown   } else {
5860ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
5870ced3a2bSJed Brown   }
5880ced3a2bSJed Brown #endif
5890ced3a2bSJed Brown   PetscFunctionReturn(0);
5900ced3a2bSJed Brown }
591e9e4536cSHong Zhang 
5928a07c6f1SJed Brown #undef __FUNCT__
5938a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
5948a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
5958a07c6f1SJed Brown {
5968a07c6f1SJed Brown   PetscErrorCode     ierr;
5978a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5988a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
5998a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6008a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6018a07c6f1SJed Brown   PetscReal          afill;
6028a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6030298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6048a07c6f1SJed Brown   PetscHeap          h;
6058a07c6f1SJed Brown   PetscBT            bt;
6068a07c6f1SJed Brown 
6078a07c6f1SJed Brown   PetscFunctionBegin;
608cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6098a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6108a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
6118a07c6f1SJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6128a07c6f1SJed Brown   ci[0] = 0;
6138a07c6f1SJed Brown 
6148a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
6158a07c6f1SJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6162205254eSKarl Rupp 
6178a07c6f1SJed Brown   current_space = free_space;
6188a07c6f1SJed Brown 
6198a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
6208a07c6f1SJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
6218a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6228a07c6f1SJed Brown 
6238a07c6f1SJed Brown   /* Determine ci and cj */
6248a07c6f1SJed Brown   for (i=0; i<am; i++) {
6258a07c6f1SJed 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 */
6268a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6278a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6288a07c6f1SJed Brown     ci[i+1] = ci[i];
6298a07c6f1SJed Brown     /* Populate the min heap */
6308a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6318a07c6f1SJed Brown       PetscInt brow = acol[j];
6328a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6338a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6348a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6358a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6368a07c6f1SJed Brown           bb[j]++;
6378a07c6f1SJed Brown           break;
6388a07c6f1SJed Brown         }
6398a07c6f1SJed Brown       }
6408a07c6f1SJed Brown     }
6418a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6428a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6438a07c6f1SJed Brown     while (j >= 0) {
6448a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6450298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
6468a07c6f1SJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
6478a07c6f1SJed Brown         ndouble++;
6488a07c6f1SJed Brown       }
6498a07c6f1SJed Brown       *(current_space->array++) = col;
6508a07c6f1SJed Brown       current_space->local_used++;
6518a07c6f1SJed Brown       current_space->local_remaining--;
6528a07c6f1SJed Brown       ci[i+1]++;
6538a07c6f1SJed Brown 
6548a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6558a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6568a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6578a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6588a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6598a07c6f1SJed Brown           bb[j]++;
6608a07c6f1SJed Brown           break;
6618a07c6f1SJed Brown         }
6628a07c6f1SJed Brown       }
6638a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6648a07c6f1SJed Brown     }
6658a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6668a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6678a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6688a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6698a07c6f1SJed Brown     }
6708a07c6f1SJed Brown   }
6718a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6728a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6738a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6748a07c6f1SJed Brown 
6758a07c6f1SJed Brown   /* Column indices are in the list of free space */
6768a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6778a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
6788a07c6f1SJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6798a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6808a07c6f1SJed Brown 
6818a07c6f1SJed Brown   /* put together the new symbolic matrix */
682ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
6832205254eSKarl Rupp 
6848a07c6f1SJed Brown   (*C)->rmap->bs = A->rmap->bs;
6858a07c6f1SJed Brown   (*C)->cmap->bs = B->cmap->bs;
6868a07c6f1SJed Brown 
6878a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6888a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6898a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
6908a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
6918a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
6928a07c6f1SJed Brown   c->nonew   = 0;
69326fbe8dcSKarl Rupp 
69489d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
6958a07c6f1SJed Brown 
6968a07c6f1SJed Brown   /* set MatInfo */
6978a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6988a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
6998a07c6f1SJed Brown   c->maxnz                     = ci[am];
7008a07c6f1SJed Brown   c->nz                        = ci[am];
7018a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7028a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7038a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7048a07c6f1SJed Brown 
7058a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7068a07c6f1SJed Brown   if (ci[am]) {
7078a07c6f1SJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
7088a07c6f1SJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
7098a07c6f1SJed Brown   } else {
7108a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7118a07c6f1SJed Brown   }
7128a07c6f1SJed Brown #endif
7138a07c6f1SJed Brown   PetscFunctionReturn(0);
7148a07c6f1SJed Brown }
7158a07c6f1SJed Brown 
716cd093f1dSJed Brown #undef __FUNCT__
71758cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
718cd093f1dSJed Brown /* concatenate unique entries and then sort */
71958cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
720cd093f1dSJed Brown {
721cd093f1dSJed Brown   PetscErrorCode     ierr;
722cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
723cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
724cd093f1dSJed Brown   PetscInt           *ci,*cj;
725cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
726cd093f1dSJed Brown   PetscReal          afill;
727cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
728cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
729cd093f1dSJed Brown   char               *seen;
730cd093f1dSJed Brown 
731cd093f1dSJed Brown   PetscFunctionBegin;
732cd093f1dSJed Brown   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
733cd093f1dSJed Brown   ci[0] = 0;
734cd093f1dSJed Brown 
735cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
736cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
737cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
738cd093f1dSJed Brown   ierr = PetscMalloc(bn*sizeof(char),&seen);CHKERRQ(ierr);
739cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
740cd093f1dSJed Brown 
741cd093f1dSJed Brown   /* Determine ci and cj */
742cd093f1dSJed Brown   for (i=0; i<am; i++) {
743cd093f1dSJed 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 */
744cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
745cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
746cd093f1dSJed Brown     /* Pack segrow */
747cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
748cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
749cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
750cd093f1dSJed Brown         PetscInt bcol = bj[k];
751cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
752cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
753cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
754cd093f1dSJed Brown           *slot = bcol;
755cd093f1dSJed Brown           seen[bcol] = 1;
756cd093f1dSJed Brown           packlen++;
757cd093f1dSJed Brown         }
758cd093f1dSJed Brown       }
759cd093f1dSJed Brown     }
760cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
761cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
762cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
763cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
764cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
765cd093f1dSJed Brown   }
766cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
767cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
768cd093f1dSJed Brown 
769cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
770cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
771cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
772cd093f1dSJed Brown 
773cd093f1dSJed Brown   /* put together the new symbolic matrix */
774cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
775cd093f1dSJed Brown 
776cd093f1dSJed Brown   (*C)->rmap->bs = A->rmap->bs;
777cd093f1dSJed Brown   (*C)->cmap->bs = B->cmap->bs;
778cd093f1dSJed Brown 
779cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
780cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
781cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
782cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
783cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
784cd093f1dSJed Brown   c->nonew   = 0;
785cd093f1dSJed Brown 
786cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
787cd093f1dSJed Brown 
788cd093f1dSJed Brown   /* set MatInfo */
789cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
790cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
791cd093f1dSJed Brown   c->maxnz                     = ci[am];
792cd093f1dSJed Brown   c->nz                        = ci[am];
793cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
794cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
795cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
796cd093f1dSJed Brown 
797cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
798cd093f1dSJed Brown   if (ci[am]) {
799cd093f1dSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
800cd093f1dSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
801cd093f1dSJed Brown   } else {
802cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
803cd093f1dSJed Brown   }
804cd093f1dSJed Brown #endif
805cd093f1dSJed Brown   PetscFunctionReturn(0);
806cd093f1dSJed Brown }
807cd093f1dSJed Brown 
808d2853540SHong Zhang /* This routine is not used. Should be removed! */
809bc011b1eSHong Zhang #undef __FUNCT__
8106fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
8116fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8125df89d91SHong Zhang {
813bc011b1eSHong Zhang   PetscErrorCode ierr;
814bc011b1eSHong Zhang 
815bc011b1eSHong Zhang   PetscFunctionBegin;
816bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
8173ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
8186fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8193ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
820bc011b1eSHong Zhang   }
8213ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
8226fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8233ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
824bc011b1eSHong Zhang   PetscFunctionReturn(0);
825bc011b1eSHong Zhang }
826bc011b1eSHong Zhang 
827bc011b1eSHong Zhang #undef __FUNCT__
828e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
829e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
8302128a86cSHong Zhang {
8312128a86cSHong Zhang   PetscErrorCode      ierr;
832e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
8332128a86cSHong Zhang 
8342128a86cSHong Zhang   PetscFunctionBegin;
835b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
8362128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
8372128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
8382128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
8392128a86cSHong Zhang   PetscFunctionReturn(0);
8402128a86cSHong Zhang }
8412128a86cSHong Zhang 
8422128a86cSHong Zhang #undef __FUNCT__
8432128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
8442128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8452128a86cSHong Zhang {
8462128a86cSHong Zhang   PetscErrorCode      ierr;
8472128a86cSHong Zhang   PetscContainer      container;
8480298fd71SBarry Smith   Mat_MatMatTransMult *multtrans=NULL;
8492128a86cSHong Zhang 
8502128a86cSHong Zhang   PetscFunctionBegin;
851e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr);
8522128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
8532128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr);
8542205254eSKarl Rupp 
8552128a86cSHong Zhang   A->ops->destroy = multtrans->destroy;
8562128a86cSHong Zhang   if (A->ops->destroy) {
8572128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
8582128a86cSHong Zhang   }
859e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
8602128a86cSHong Zhang   PetscFunctionReturn(0);
8612128a86cSHong Zhang }
8622128a86cSHong Zhang 
8632128a86cSHong Zhang #undef __FUNCT__
8646fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
8656fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
866bc011b1eSHong Zhang {
8675df89d91SHong Zhang   PetscErrorCode      ierr;
86881d82fe4SHong Zhang   Mat                 Bt;
86981d82fe4SHong Zhang   PetscInt            *bti,*btj;
870e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
8712128a86cSHong Zhang   PetscContainer      container;
872d2853540SHong Zhang 
87381d82fe4SHong Zhang   PetscFunctionBegin;
87481d82fe4SHong Zhang   /* create symbolic Bt */
87581d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8760298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
8772205254eSKarl Rupp 
878c3e5699bSSatish Balay   Bt->rmap->bs = A->cmap->bs;
879c3e5699bSSatish Balay   Bt->cmap->bs = B->cmap->bs;
88081d82fe4SHong Zhang 
88181d82fe4SHong Zhang   /* get symbolic C=A*Bt */
88281d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
88381d82fe4SHong Zhang 
8842128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
885e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
8862128a86cSHong Zhang 
8872128a86cSHong Zhang   /* attach the supporting struct to C */
8882128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
8892128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
890e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
891e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
8922128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
8932128a86cSHong Zhang 
8942128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
8952128a86cSHong Zhang   multtrans->destroy     = (*C)->ops->destroy;
8962128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8972128a86cSHong Zhang 
8980298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-matmattransmult_color",&multtrans->usecoloring,NULL);CHKERRQ(ierr);
8992128a86cSHong Zhang   if (multtrans->usecoloring) {
900b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
901b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
9022128a86cSHong Zhang     ISColoring           iscoloring;
9032128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
904e8354b3eSHong Zhang 
9054614e006SHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
906b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
9072205254eSKarl Rupp 
9082128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
9092205254eSKarl Rupp 
9102128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
9112128a86cSHong Zhang 
9122128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
9132128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
9142128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9152128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
9160298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9172205254eSKarl Rupp 
9182128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
9192128a86cSHong Zhang     multtrans->Bt_den   = Bt_dense;
9202128a86cSHong Zhang 
9212128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9222128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9232128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9240298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9252205254eSKarl Rupp 
9262128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
9272128a86cSHong Zhang     multtrans->ABt_den  = C_dense;
928f94ccd6cSHong Zhang 
929f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9301ce71dffSSatish Balay     {
931f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
93222d28d08SBarry 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);
9331ce71dffSSatish Balay     }
934f94ccd6cSHong Zhang #endif
9352128a86cSHong Zhang   }
936e99be685SHong Zhang   /* clean up */
937e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
938e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9395df89d91SHong Zhang   PetscFunctionReturn(0);
9405df89d91SHong Zhang }
9415df89d91SHong Zhang 
9426973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
9435df89d91SHong Zhang #undef __FUNCT__
9446fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9456fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9465df89d91SHong Zhang {
9475df89d91SHong Zhang   PetscErrorCode      ierr;
9485df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
949e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
95081d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9515df89d91SHong Zhang   PetscLogDouble      flops=0.0;
952aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
953e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
9542128a86cSHong Zhang   PetscContainer      container;
9556973769bSHong Zhang #if defined(USE_ARRAY)
9566973769bSHong Zhang   MatScalar *spdot;
9576973769bSHong Zhang #endif
9585df89d91SHong Zhang 
9595df89d91SHong Zhang   PetscFunctionBegin;
96058ed3ceaSHong Zhang   /* clear old values in C */
96158ed3ceaSHong Zhang   if (!c->a) {
96258ed3ceaSHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
96358ed3ceaSHong Zhang     c->a      = ca;
96458ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
96558ed3ceaSHong Zhang   } else {
96658ed3ceaSHong Zhang     ca =  c->a;
96758ed3ceaSHong Zhang   }
96858ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
96958ed3ceaSHong Zhang 
970e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr);
9712128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
9722128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr);
9732128a86cSHong Zhang   if (multtrans->usecoloring) {
974b9af6bddSHong Zhang     MatTransposeColoring matcoloring = multtrans->matcoloring;
975c8db22f6SHong Zhang     Mat                  Bt_dense;
976c8db22f6SHong Zhang     PetscInt             m,n;
977a0b95ffbSSatish Balay     Mat                  C_dense = multtrans->ABt_den;
978*b3b4fc5aSHong Zhang     PetscLogDouble       t0,t1,t2,t3,Bt_den,C_den,C_sp;
979c8db22f6SHong Zhang 
980b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
981*b3b4fc5aSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
982*b3b4fc5aSHong Zhang     Bt_dense = multtrans->Bt_den;
983b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
984*b3b4fc5aSHong Zhang     ierr = PetscGetTime(&t1);CHKERRQ(ierr);
985*b3b4fc5aSHong Zhang     Bt_den = t1 - t0;
986c8db22f6SHong Zhang 
987c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9882128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
989*b3b4fc5aSHong Zhang     ierr = PetscGetTime(&t2);CHKERRQ(ierr);
990*b3b4fc5aSHong Zhang     C_den = t2 - t1;
991c8db22f6SHong Zhang 
9922128a86cSHong Zhang     /* Recover C from C_dense */
993b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
994*b3b4fc5aSHong Zhang     ierr = PetscGetTime(&t3);CHKERRQ(ierr);
995*b3b4fc5aSHong Zhang     C_sp = t3 - t2;
996*b3b4fc5aSHong Zhang #if defined(PETSC_USE_INFO)
997*b3b4fc5aSHong 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);
998*b3b4fc5aSHong Zhang     ierr     = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
999*b3b4fc5aSHong Zhang     ierr = PetscInfo4(C,"Bt_den: %d x %d, B: %d x %d\n",m,n,m,C->cmap->n);CHKERRQ(ierr);
1000*b3b4fc5aSHong Zhang #endif
1001c8db22f6SHong Zhang     PetscFunctionReturn(0);
1002c8db22f6SHong Zhang   }
1003c8db22f6SHong Zhang 
10046973769bSHong Zhang #if defined(USE_ARRAY)
10056973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
1006e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
1007e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
10086973769bSHong Zhang #endif
10096973769bSHong Zhang 
10101fa1209cSHong Zhang   for (i=0; i<cm; i++) {
101181d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
101281d82fe4SHong Zhang     acol = aj + ai[i];
10136973769bSHong Zhang     aval = aa + ai[i];
10141fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
10151fa1209cSHong Zhang     ccol = cj + ci[i];
10166973769bSHong Zhang     cval = ca + ci[i];
10171fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
101881d82fe4SHong Zhang       brow = ccol[j];
101981d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
102081d82fe4SHong Zhang       bcol = bj + bi[brow];
10216973769bSHong Zhang       bval = ba + bi[brow];
10226973769bSHong Zhang 
102381d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
10246973769bSHong Zhang #if defined(USE_ARRAY)
10256973769bSHong Zhang       /* put ba to spdot array */
10266973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
10276973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
10286973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++) {
10296973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
10306973769bSHong Zhang       }
10316973769bSHong Zhang       /* zero spdot array */
10326973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
10336973769bSHong Zhang #else
103481d82fe4SHong Zhang       nexta = 0; nextb = 0;
103581d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
10367b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
103781d82fe4SHong Zhang         if (nexta == anzi) break;
10387b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
103981d82fe4SHong Zhang         if (nextb == bnzj) break;
104081d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
10416973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
104281d82fe4SHong Zhang           nexta++; nextb++;
104381d82fe4SHong Zhang           flops += 2;
10441fa1209cSHong Zhang         }
10451fa1209cSHong Zhang       }
10466973769bSHong Zhang #endif
104781d82fe4SHong Zhang     }
104881d82fe4SHong Zhang   }
104981d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105081d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105181d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10526973769bSHong Zhang #if defined(USE_ARRAY)
10536973769bSHong Zhang   ierr = PetscFree(spdot);
10546973769bSHong Zhang #endif
10555df89d91SHong Zhang   PetscFunctionReturn(0);
10565df89d91SHong Zhang }
10575df89d91SHong Zhang 
10585df89d91SHong Zhang #undef __FUNCT__
105975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
10600adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10610adebc6cSBarry Smith {
10625df89d91SHong Zhang   PetscErrorCode ierr;
10635df89d91SHong Zhang 
10645df89d91SHong Zhang   PetscFunctionBegin;
10655df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
106675648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
10675df89d91SHong Zhang   }
106875648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
10695df89d91SHong Zhang   PetscFunctionReturn(0);
10705df89d91SHong Zhang }
10715df89d91SHong Zhang 
10725df89d91SHong Zhang #undef __FUNCT__
107375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
107475648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10755df89d91SHong Zhang {
1076bc011b1eSHong Zhang   PetscErrorCode ierr;
1077bc011b1eSHong Zhang   Mat            At;
107838baddfdSBarry Smith   PetscInt       *ati,*atj;
1079bc011b1eSHong Zhang 
1080bc011b1eSHong Zhang   PetscFunctionBegin;
1081bc011b1eSHong Zhang   /* create symbolic At */
1082bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10830298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
10842205254eSKarl Rupp 
1085a2f3521dSMark F. Adams   At->rmap->bs = A->cmap->bs;
1086a2f3521dSMark F. Adams   At->cmap->bs = B->cmap->bs;
1087bc011b1eSHong Zhang 
1088bc011b1eSHong Zhang   /* get symbolic C=At*B */
1089bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1090bc011b1eSHong Zhang 
1091bc011b1eSHong Zhang   /* clean up */
10926bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1093bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1094bc011b1eSHong Zhang   PetscFunctionReturn(0);
1095bc011b1eSHong Zhang }
1096bc011b1eSHong Zhang 
1097bc011b1eSHong Zhang #undef __FUNCT__
109875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
109975648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1100bc011b1eSHong Zhang {
1101bc011b1eSHong Zhang   PetscErrorCode ierr;
11020fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1103d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1104d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1105d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1106aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1107bc011b1eSHong Zhang 
1108bc011b1eSHong Zhang   PetscFunctionBegin;
1109aa1aec99SHong Zhang   if (!c->a) {
1110aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
11112205254eSKarl Rupp 
1112aa1aec99SHong Zhang     c->a      = ca;
1113aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1114aa1aec99SHong Zhang   } else {
1115aa1aec99SHong Zhang     ca = c->a;
1116aa1aec99SHong Zhang   }
1117bc011b1eSHong Zhang   /* clear old values in C */
1118bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1119bc011b1eSHong Zhang 
1120bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1121bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1122bc011b1eSHong Zhang     bj   = b->j + bi[i];
1123bc011b1eSHong Zhang     ba   = b->a + bi[i];
1124bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1125bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1126bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1127bc011b1eSHong Zhang       nextb = 0;
11280fbc74f4SHong Zhang       crow  = *aj++;
1129bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1130bc011b1eSHong Zhang       caj   = ca + ci[crow];
1131bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1132bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
11330fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
11340fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1135bc011b1eSHong Zhang           nextb++;
1136bc011b1eSHong Zhang         }
1137bc011b1eSHong Zhang       }
1138bc011b1eSHong Zhang       flops += 2*bnzi;
11390fbc74f4SHong Zhang       aa++;
1140bc011b1eSHong Zhang     }
1141bc011b1eSHong Zhang   }
1142bc011b1eSHong Zhang 
1143bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1144bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1145bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1146bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1147bc011b1eSHong Zhang   PetscFunctionReturn(0);
1148bc011b1eSHong Zhang }
11499123193aSHong Zhang 
11509123193aSHong Zhang #undef __FUNCT__
11519123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
11529123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11539123193aSHong Zhang {
11549123193aSHong Zhang   PetscErrorCode ierr;
11559123193aSHong Zhang 
11569123193aSHong Zhang   PetscFunctionBegin;
11579123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
11583ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11599123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11603ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11619123193aSHong Zhang   }
11623ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11639123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11644614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11659123193aSHong Zhang   PetscFunctionReturn(0);
11669123193aSHong Zhang }
1167edd81eecSMatthew Knepley 
11689123193aSHong Zhang #undef __FUNCT__
11699123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11709123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11719123193aSHong Zhang {
11729123193aSHong Zhang   PetscErrorCode ierr;
11739123193aSHong Zhang 
11749123193aSHong Zhang   PetscFunctionBegin;
11755a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11762205254eSKarl Rupp 
1177d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11789123193aSHong Zhang   PetscFunctionReturn(0);
11799123193aSHong Zhang }
11809123193aSHong Zhang 
11819123193aSHong Zhang #undef __FUNCT__
11829123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11839123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11849123193aSHong Zhang {
1185f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11869123193aSHong Zhang   PetscErrorCode ierr;
1187dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1188dd6ea824SBarry Smith   MatScalar      *aa;
1189d0f46423SBarry Smith   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
1190f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
11919123193aSHong Zhang 
11929123193aSHong Zhang   PetscFunctionBegin;
1193f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1194e32f2f54SBarry 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);
1195e32f2f54SBarry 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);
1196e32f2f54SBarry 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);
11978c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
11988c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1199f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1200f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1201f32d5d43SBarry Smith     colam = col*am;
1202f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1203f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1204f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1205f32d5d43SBarry Smith       aj = a->j + a->i[i];
1206f32d5d43SBarry Smith       aa = a->a + a->i[i];
1207f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1208f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
1209f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
1210f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
1211f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
12129123193aSHong Zhang       }
1213f32d5d43SBarry Smith       c[colam + i]       = r1;
1214f32d5d43SBarry Smith       c[colam + am + i]  = r2;
1215f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
1216f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
1217f32d5d43SBarry Smith     }
1218f32d5d43SBarry Smith     b1 += bm4;
1219f32d5d43SBarry Smith     b2 += bm4;
1220f32d5d43SBarry Smith     b3 += bm4;
1221f32d5d43SBarry Smith     b4 += bm4;
1222f32d5d43SBarry Smith   }
1223f32d5d43SBarry Smith   for (; col<cn; col++) {     /* over extra columns of C */
1224f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1225f32d5d43SBarry Smith       r1 = 0.0;
1226f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1227f32d5d43SBarry Smith       aj = a->j + a->i[i];
1228f32d5d43SBarry Smith       aa = a->a + a->i[i];
1229f32d5d43SBarry Smith 
1230f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1231f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
1232f32d5d43SBarry Smith       }
1233f32d5d43SBarry Smith       c[col*am + i] = r1;
1234f32d5d43SBarry Smith     }
1235f32d5d43SBarry Smith     b1 += bm;
1236f32d5d43SBarry Smith   }
1237dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
12388c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
12398c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1240f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1241f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1242f32d5d43SBarry Smith   PetscFunctionReturn(0);
1243f32d5d43SBarry Smith }
1244f32d5d43SBarry Smith 
1245f32d5d43SBarry Smith /*
12464324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1247f32d5d43SBarry Smith */
1248f32d5d43SBarry Smith #undef __FUNCT__
1249f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1250f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1251f32d5d43SBarry Smith {
1252f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1253f32d5d43SBarry Smith   PetscErrorCode ierr;
1254dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1255dd6ea824SBarry Smith   MatScalar      *aa;
1256d0f46423SBarry 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;
12574324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1258f32d5d43SBarry Smith 
1259f32d5d43SBarry Smith   PetscFunctionBegin;
1260f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
12618c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12628c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1263f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12644324174eSBarry Smith 
12654324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12664324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12674324174eSBarry Smith       colam = col*am;
12684324174eSBarry Smith       arm   = a->compressedrow.nrows;
12694324174eSBarry Smith       ii    = a->compressedrow.i;
12704324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12714324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12724324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12734324174eSBarry Smith         n  = ii[i+1] - ii[i];
12744324174eSBarry Smith         aj = a->j + ii[i];
12754324174eSBarry Smith         aa = a->a + ii[i];
12764324174eSBarry Smith         for (j=0; j<n; j++) {
12774324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12784324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12794324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12804324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12814324174eSBarry Smith         }
12824324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12834324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12844324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12854324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12864324174eSBarry Smith       }
12874324174eSBarry Smith       b1 += bm4;
12884324174eSBarry Smith       b2 += bm4;
12894324174eSBarry Smith       b3 += bm4;
12904324174eSBarry Smith       b4 += bm4;
12914324174eSBarry Smith     }
12924324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12934324174eSBarry Smith       colam = col*am;
12944324174eSBarry Smith       arm   = a->compressedrow.nrows;
12954324174eSBarry Smith       ii    = a->compressedrow.i;
12964324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12974324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12984324174eSBarry Smith         r1 = 0.0;
12994324174eSBarry Smith         n  = ii[i+1] - ii[i];
13004324174eSBarry Smith         aj = a->j + ii[i];
13014324174eSBarry Smith         aa = a->a + ii[i];
13024324174eSBarry Smith 
13034324174eSBarry Smith         for (j=0; j<n; j++) {
13044324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
13054324174eSBarry Smith         }
1306a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
13074324174eSBarry Smith       }
13084324174eSBarry Smith       b1 += bm;
13094324174eSBarry Smith     }
13104324174eSBarry Smith   } else {
1311f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1312f32d5d43SBarry Smith       colam = col*am;
1313f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1314f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1315f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1316f32d5d43SBarry Smith         aj = a->j + a->i[i];
1317f32d5d43SBarry Smith         aa = a->a + a->i[i];
1318f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1319f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1320f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1321f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1322f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1323f32d5d43SBarry Smith         }
1324f32d5d43SBarry Smith         c[colam + i]       += r1;
1325f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1326f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1327f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1328f32d5d43SBarry Smith       }
1329f32d5d43SBarry Smith       b1 += bm4;
1330f32d5d43SBarry Smith       b2 += bm4;
1331f32d5d43SBarry Smith       b3 += bm4;
1332f32d5d43SBarry Smith       b4 += bm4;
1333f32d5d43SBarry Smith     }
1334f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1335a2ea699eSBarry Smith       colam = col*am;
1336f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1337f32d5d43SBarry Smith         r1 = 0.0;
1338f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1339f32d5d43SBarry Smith         aj = a->j + a->i[i];
1340f32d5d43SBarry Smith         aa = a->a + a->i[i];
1341f32d5d43SBarry Smith 
1342f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1343f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1344f32d5d43SBarry Smith         }
1345a2ea699eSBarry Smith         c[colam + i] += r1;
1346f32d5d43SBarry Smith       }
1347f32d5d43SBarry Smith       b1 += bm;
1348f32d5d43SBarry Smith     }
13494324174eSBarry Smith   }
1350dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13518c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
13528c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13539123193aSHong Zhang   PetscFunctionReturn(0);
13549123193aSHong Zhang }
1355b1683b59SHong Zhang 
1356b1683b59SHong Zhang #undef __FUNCT__
1357b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1358b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1359c8db22f6SHong Zhang {
1360c8db22f6SHong Zhang   PetscErrorCode ierr;
13612f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13622f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13632f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13642f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13652f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13662f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1367c8db22f6SHong Zhang 
1368c8db22f6SHong Zhang   PetscFunctionBegin;
13692f5f1f90SHong Zhang   btval_den=btdense->v;
13702f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13712f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13722f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13732f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1374d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13752f5f1f90SHong Zhang       btcol = bj + bi[col];
13762f5f1f90SHong Zhang       btval = ba + bi[col];
13772f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1378d2853540SHong Zhang       for (j=0; j<anz; j++) {
13792f5f1f90SHong Zhang         brow            = btcol[j];
13802f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1381c8db22f6SHong Zhang       }
1382c8db22f6SHong Zhang     }
13832f5f1f90SHong Zhang     btval_den += m;
1384c8db22f6SHong Zhang   }
1385c8db22f6SHong Zhang   PetscFunctionReturn(0);
1386c8db22f6SHong Zhang }
1387c8db22f6SHong Zhang 
1388c8db22f6SHong Zhang #undef __FUNCT__
1389b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1390b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13918972f759SHong Zhang {
13928972f759SHong Zhang   PetscErrorCode ierr;
1393b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
13942f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1395b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
13962f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
13978972f759SHong Zhang 
13988972f759SHong Zhang   PetscFunctionBegin;
13990298fd71SBarry Smith   ierr   = MatGetLocalSize(Csp,&m,NULL);CHKERRQ(ierr);
1400a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1401b2d2b04fSHong Zhang   cp_den = ca_den;
14022f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
14032f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
14042f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
14052f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
14062f5f1f90SHong Zhang     for (l=0; l<nrows; l++) {
14072f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1408b2d2b04fSHong Zhang     }
1409b2d2b04fSHong Zhang     cp_den += m;
1410b2d2b04fSHong Zhang   }
1411a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
14128972f759SHong Zhang   PetscFunctionReturn(0);
14138972f759SHong Zhang }
14148972f759SHong Zhang 
1415e99be685SHong Zhang /*
1416e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1417e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1418e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1419e99be685SHong Zhang  */
1420e99be685SHong Zhang #undef __FUNCT__
1421e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
14221a83f524SJed 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)
1423e99be685SHong Zhang {
1424e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1425e99be685SHong Zhang   PetscErrorCode ierr;
1426e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1427e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1428e99be685SHong Zhang   PetscInt       *cspidx;
1429e99be685SHong Zhang 
1430e99be685SHong Zhang   PetscFunctionBegin;
1431e99be685SHong Zhang   *nn = n;
1432e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1433e99be685SHong Zhang   if (symmetric) {
1434ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
14351a83f524SJed Brown     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr);
1436e99be685SHong Zhang   } else {
1437e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1438e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1439e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1440e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1441e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1442e99be685SHong Zhang     jj   = a->j;
1443e99be685SHong Zhang     for (i=0; i<nz; i++) {
1444e99be685SHong Zhang       collengths[jj[i]]++;
1445e99be685SHong Zhang     }
1446e99be685SHong Zhang     cia[0] = oshift;
1447e99be685SHong Zhang     for (i=0; i<n; i++) {
1448e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1449e99be685SHong Zhang     }
1450e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1451e99be685SHong Zhang     jj   = a->j;
1452e99be685SHong Zhang     for (row=0; row<m; row++) {
1453e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1454e99be685SHong Zhang       for (i=0; i<mr; i++) {
1455e99be685SHong Zhang         col = *jj++;
14562205254eSKarl Rupp 
1457e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1458e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift]  = row + oshift;
1459e99be685SHong Zhang       }
1460e99be685SHong Zhang     }
1461e99be685SHong Zhang     ierr   = PetscFree(collengths);CHKERRQ(ierr);
1462e99be685SHong Zhang     *ia    = cia; *ja = cja;
1463e99be685SHong Zhang     *spidx = cspidx;
1464e99be685SHong Zhang   }
1465e99be685SHong Zhang   PetscFunctionReturn(0);
1466e99be685SHong Zhang }
1467e99be685SHong Zhang 
1468e99be685SHong Zhang #undef __FUNCT__
1469e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
14701a83f524SJed 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)
1471e99be685SHong Zhang {
1472e99be685SHong Zhang   PetscErrorCode ierr;
1473e99be685SHong Zhang 
1474e99be685SHong Zhang   PetscFunctionBegin;
1475e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1476e99be685SHong Zhang 
1477e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1478e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1479e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1480e99be685SHong Zhang   PetscFunctionReturn(0);
1481e99be685SHong Zhang }
1482e99be685SHong Zhang 
14838972f759SHong Zhang #undef __FUNCT__
1484b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1485b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1486b1683b59SHong Zhang {
1487b1683b59SHong Zhang   PetscErrorCode ierr;
14881a83f524SJed Brown   PetscInt       i,n,nrows,N,j,k,m,ncols,col,cm;
14891a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1490b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1491b1683b59SHong Zhang   IS             *isa;
1492b1683b59SHong Zhang   PetscBool      flg1,flg2;
1493b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1494e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1495d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1496b1683b59SHong Zhang 
1497b1683b59SHong Zhang   PetscFunctionBegin;
1498b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1499e99be685SHong Zhang 
1500b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1501251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1502251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1503b1683b59SHong Zhang   if (flg1 || flg2) {
1504b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1505b1683b59SHong Zhang   }
1506b1683b59SHong Zhang 
1507b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1508b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1509b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1510b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1511b1683b59SHong Zhang   c->rstart = 0;
1512b1683b59SHong Zhang 
1513b1683b59SHong Zhang   c->ncolors = nis;
1514b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1515b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1516d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1517d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
15182205254eSKarl Rupp 
1519d2853540SHong Zhang   colorforrow[0]    = 0;
1520d2853540SHong Zhang   rows_i            = rows;
1521d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1522b1683b59SHong Zhang 
1523d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1524d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
15252205254eSKarl Rupp 
1526d2853540SHong Zhang   colorforcol[0] = 0;
1527d2853540SHong Zhang   columns_i      = columns;
1528d2853540SHong Zhang 
15290298fd71SBarry Smith   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); /* column-wise storage of mat */
1530b1683b59SHong Zhang 
1531b28d80bdSHong Zhang   cm   = c->m;
1532b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1533e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1534b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1535b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1536b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
15372205254eSKarl Rupp 
1538b1683b59SHong Zhang     c->ncolumns[i] = n;
1539b1683b59SHong Zhang     if (n) {
1540d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1541b1683b59SHong Zhang     }
1542d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1543d2853540SHong Zhang     columns_i       += n;
1544b1683b59SHong Zhang 
1545b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1546e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1547e99be685SHong Zhang 
1548b1683b59SHong Zhang     /* loop over columns*/
1549b1683b59SHong Zhang     for (j=0; j<n; j++) {
1550b1683b59SHong Zhang       col     = is[j];
1551d2853540SHong Zhang       row_idx = cj + ci[col];
1552b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1553b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1554b1683b59SHong Zhang       for (k=0; k<m; k++) {
1555e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1556d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1557b1683b59SHong Zhang       }
1558b1683b59SHong Zhang     }
1559b1683b59SHong Zhang     /* count the number of hits */
1560b1683b59SHong Zhang     nrows = 0;
1561e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1562b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1563b1683b59SHong Zhang     }
1564b1683b59SHong Zhang     c->nrows[i]      = nrows;
1565d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1566d2853540SHong Zhang 
1567b1683b59SHong Zhang     nrows = 0;
1568e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1569b1683b59SHong Zhang       if (rowhit[j]) {
1570d2853540SHong Zhang         rows_i[nrows]            = j;
1571e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1572b1683b59SHong Zhang         nrows++;
1573b1683b59SHong Zhang       }
1574b1683b59SHong Zhang     }
1575b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1576d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1577b1683b59SHong Zhang   }
15780298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1579b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1580b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1581d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1582d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1583d2853540SHong Zhang #endif
1584b28d80bdSHong Zhang 
1585d2853540SHong Zhang   c->colorforrow     = colorforrow;
1586d2853540SHong Zhang   c->rows            = rows;
1587d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1588d2853540SHong Zhang   c->colorforcol     = colorforcol;
1589d2853540SHong Zhang   c->columns         = columns;
15902205254eSKarl Rupp 
1591f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1592b1683b59SHong Zhang   PetscFunctionReturn(0);
1593b1683b59SHong Zhang }
1594