xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 335efc43fa2eb4d5736843c09c28f9f3dff66c92)
1be1d678aSKris Buschelman 
2d50806bdSBarry Smith /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4d50806bdSBarry Smith           C = A * B
5d50806bdSBarry Smith */
6d50806bdSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
90ced3a2bSJed Brown #include <../src/mat/utils/petscheap.h>
10c6db04a5SJed Brown #include <petscbt.h>
1107475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h>
1226be7272SHong Zhang #include <petsctime.h>
137bab7c10SHong Zhang 
1458cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
15cd093f1dSJed Brown 
166284ec50SHong Zhang #undef __FUNCT__
175c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1838baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1938baddfdSBarry Smith {
20dfbe8321SBarry Smith   PetscErrorCode ierr;
216540a9cdSHong Zhang   const char     *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
226540a9cdSHong Zhang   PetscInt       alg=0; /* set default algorithm */
235c66b693SKris Buschelman 
245c66b693SKris Buschelman   PetscFunctionBegin;
2526be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
26152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
276540a9cdSHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,6,algTypes[0],&alg,NULL);CHKERRQ(ierr);
28d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
293ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
306540a9cdSHong Zhang     switch (alg) {
316540a9cdSHong Zhang     case 1:
32aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
336540a9cdSHong Zhang       break;
346540a9cdSHong Zhang     case 2:
356540a9cdSHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
366540a9cdSHong Zhang       break;
376540a9cdSHong Zhang     case 3:
380ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
396540a9cdSHong Zhang       break;
406540a9cdSHong Zhang     case 4:
418a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
426540a9cdSHong Zhang       break;
436540a9cdSHong Zhang     case 5:
4458cf0668SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
456540a9cdSHong Zhang       break;
466540a9cdSHong Zhang     default:
4726be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
486540a9cdSHong Zhang      break;
4925023028SHong Zhang     }
503ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
5126be0446SHong Zhang   }
525c913ed7SHong Zhang 
533ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
5401e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
553ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
565c66b693SKris Buschelman   PetscFunctionReturn(0);
575c66b693SKris Buschelman }
581c24bd37SHong Zhang 
59e9e4536cSHong Zhang #undef __FUNCT__
6058cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed"
6158cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
62b561aa9dSHong Zhang {
63b561aa9dSHong Zhang   PetscErrorCode     ierr;
64b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
65971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
66c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
67b561aa9dSHong Zhang   PetscReal          afill;
68fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
69fb908031SHong Zhang   PetscBT            lnkbt;
700298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
71b561aa9dSHong Zhang 
72b561aa9dSHong Zhang   PetscFunctionBegin;
73bd958071SHong Zhang   /* Get ci and cj */
74bd958071SHong Zhang   /*---------------*/
75fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
76fb908031SHong Zhang   /* free space for accumulating nonzero column info */
77fb908031SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
78fb908031SHong Zhang   ci[0] = 0;
79fb908031SHong Zhang 
80fb908031SHong Zhang   /* create and initialize a linked list */
81fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
82fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
83fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
84fb908031SHong Zhang 
85fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
86fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
872205254eSKarl Rupp 
88fb908031SHong Zhang   current_space = free_space;
89fb908031SHong Zhang 
90fb908031SHong Zhang   /* Determine ci and cj */
91fb908031SHong Zhang   for (i=0; i<am; i++) {
92fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
93fb908031SHong Zhang     aj   = a->j + ai[i];
94fb908031SHong Zhang     for (j=0; j<anzi; j++) {
95fb908031SHong Zhang       brow = aj[j];
96fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
97fb908031SHong Zhang       bj   = b->j + bi[brow];
98fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
99fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
100fb908031SHong Zhang     }
101fb908031SHong Zhang     cnzi = lnk[0];
102fb908031SHong Zhang 
103fb908031SHong Zhang     /* If free space is not available, make more free space */
104fb908031SHong Zhang     /* Double the amount of total space in the list */
105fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
106fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
107fb908031SHong Zhang       ndouble++;
108fb908031SHong Zhang     }
109fb908031SHong Zhang 
110fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
111fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1122205254eSKarl Rupp 
113fb908031SHong Zhang     current_space->array           += cnzi;
114fb908031SHong Zhang     current_space->local_used      += cnzi;
115fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1162205254eSKarl Rupp 
117fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
118fb908031SHong Zhang   }
119fb908031SHong Zhang 
120fb908031SHong Zhang   /* Column indices are in the list of free space */
121fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
122fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
123fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
124fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
125fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
126b561aa9dSHong Zhang 
12726be0446SHong Zhang   /* put together the new symbolic matrix */
128ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1292205254eSKarl Rupp 
130a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
131a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
13258c24d83SHong Zhang 
13358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
13458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
13558c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
136aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
137e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
13858c24d83SHong Zhang   c->nonew                  = 0;
139002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1400b7e3e3dSHong Zhang 
1417212da7cSHong Zhang   /* set MatInfo */
1427212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
143f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1447212da7cSHong Zhang   c->maxnz                     = ci[am];
1457212da7cSHong Zhang   c->nz                        = ci[am];
146fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1477212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1487212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1497212da7cSHong Zhang 
1507212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1517212da7cSHong Zhang   if (ci[am]) {
152fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
153f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
154f2b054eeSHong Zhang   } else {
155f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
156be0fcf8dSHong Zhang   }
157f2b054eeSHong Zhang #endif
15858c24d83SHong Zhang   PetscFunctionReturn(0);
15958c24d83SHong Zhang }
160d50806bdSBarry Smith 
16126be0446SHong Zhang #undef __FUNCT__
16226be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
163dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
164d50806bdSBarry Smith {
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
167d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
168d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
169d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
17038baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
171c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
172519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
173aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
174aa1aec99SHong Zhang   PetscScalar    *ab_dense;
175d50806bdSBarry Smith 
176d50806bdSBarry Smith   PetscFunctionBegin;
177aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
178aa1aec99SHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
179aa1aec99SHong Zhang     c->a      = ca;
180aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
181aa1aec99SHong Zhang   } else {
182aa1aec99SHong Zhang     ca        = c->a;
183d90d86d0SMatthew G. Knepley   }
184d90d86d0SMatthew G. Knepley   if (!c->matmult_abdense) {
185d90d86d0SMatthew G. Knepley     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
186d90d86d0SMatthew G. Knepley     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
187d90d86d0SMatthew G. Knepley     c->matmult_abdense = ab_dense;
188d90d86d0SMatthew G. Knepley   } else {
189aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
190aa1aec99SHong Zhang   }
191aa1aec99SHong Zhang 
192c124e916SHong Zhang   /* clean old values in C */
193c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
194d50806bdSBarry Smith   /* Traverse A row-wise. */
195d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
196d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
197d50806bdSBarry Smith   for (i=0; i<am; i++) {
198d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
199d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
200519eb980SHong Zhang       brow = aj[j];
201d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
202d50806bdSBarry Smith       bjj  = bj + bi[brow];
203d50806bdSBarry Smith       baj  = ba + bi[brow];
204519eb980SHong Zhang       /* perform dense axpy */
20536ec6d2dSHong Zhang       valtmp = aa[j];
206519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
20736ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
208519eb980SHong Zhang       }
209519eb980SHong Zhang       flops += 2*bnzi;
210519eb980SHong Zhang     }
211c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
212c58ca2e3SHong Zhang 
213c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
214519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
215519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
216519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
217519eb980SHong Zhang     }
218519eb980SHong Zhang     flops += cnzi;
219519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
220519eb980SHong Zhang   }
221c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
222c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
224c58ca2e3SHong Zhang   PetscFunctionReturn(0);
225c58ca2e3SHong Zhang }
226c58ca2e3SHong Zhang 
227c58ca2e3SHong Zhang #undef __FUNCT__
22825023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
22925023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
230c58ca2e3SHong Zhang {
231c58ca2e3SHong Zhang   PetscErrorCode ierr;
232c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
233c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
234c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
235c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
236c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
237c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
238c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
23936ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
240c58ca2e3SHong Zhang   PetscInt       nextb;
241c58ca2e3SHong Zhang 
242c58ca2e3SHong Zhang   PetscFunctionBegin;
243c58ca2e3SHong Zhang   /* clean old values in C */
244c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
245c58ca2e3SHong Zhang   /* Traverse A row-wise. */
246c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
247c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
248519eb980SHong Zhang   for (i=0; i<am; i++) {
249519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
250519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
251519eb980SHong Zhang     for (j=0; j<anzi; j++) {
252519eb980SHong Zhang       brow = aj[j];
253519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
254519eb980SHong Zhang       bjj  = bj + bi[brow];
255519eb980SHong Zhang       baj  = ba + bi[brow];
256519eb980SHong Zhang       /* perform sparse axpy */
25736ec6d2dSHong Zhang       valtmp = aa[j];
258c124e916SHong Zhang       nextb  = 0;
259c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
260c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
26136ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
262c124e916SHong Zhang         }
263d50806bdSBarry Smith       }
264d50806bdSBarry Smith       flops += 2*bnzi;
265d50806bdSBarry Smith     }
266519eb980SHong Zhang     aj += anzi; aa += anzi;
267519eb980SHong Zhang     cj += cnzi; ca += cnzi;
268d50806bdSBarry Smith   }
269c58ca2e3SHong Zhang 
270716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
271716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
272d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
273d50806bdSBarry Smith   PetscFunctionReturn(0);
274d50806bdSBarry Smith }
275bc011b1eSHong Zhang 
276e9e4536cSHong Zhang #undef __FUNCT__
2773c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2783c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
27925296bd5SBarry Smith {
28025296bd5SBarry Smith   PetscErrorCode     ierr;
28125296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
28225296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2833c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
28425296bd5SBarry Smith   MatScalar          *ca;
28525296bd5SBarry Smith   PetscReal          afill;
2863c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2870298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
28825296bd5SBarry Smith 
28925296bd5SBarry Smith   PetscFunctionBegin;
2903c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2913c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2923c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2933c50cad2SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2943c50cad2SHong Zhang   ci[0] = 0;
2953c50cad2SHong Zhang 
2963c50cad2SHong Zhang   /* create and initialize a linked list */
2973c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2983c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2993c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
3003c50cad2SHong Zhang 
3013c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
3023c50cad2SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
3033c50cad2SHong Zhang   current_space = free_space;
3043c50cad2SHong Zhang 
3053c50cad2SHong Zhang   /* Determine ci and cj */
3063c50cad2SHong Zhang   for (i=0; i<am; i++) {
3073c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3083c50cad2SHong Zhang     aj   = a->j + ai[i];
3093c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3103c50cad2SHong Zhang       brow = aj[j];
3113c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3123c50cad2SHong Zhang       bj   = b->j + bi[brow];
3133c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3143c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3153c50cad2SHong Zhang     }
3163c50cad2SHong Zhang     cnzi = lnk[1];
3173c50cad2SHong Zhang 
3183c50cad2SHong Zhang     /* If free space is not available, make more free space */
3193c50cad2SHong Zhang     /* Double the amount of total space in the list */
3203c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3213c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3223c50cad2SHong Zhang       ndouble++;
3233c50cad2SHong Zhang     }
3243c50cad2SHong Zhang 
3253c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3263c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3272205254eSKarl Rupp 
3283c50cad2SHong Zhang     current_space->array           += cnzi;
3293c50cad2SHong Zhang     current_space->local_used      += cnzi;
3303c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3312205254eSKarl Rupp 
3323c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3333c50cad2SHong Zhang   }
3343c50cad2SHong Zhang 
3353c50cad2SHong Zhang   /* Column indices are in the list of free space */
3363c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3373c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3383c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3393c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3403c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
34125296bd5SBarry Smith 
34225296bd5SBarry Smith   /* Allocate space for ca */
34325296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
34425296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
34525296bd5SBarry Smith 
34625296bd5SBarry Smith   /* put together the new symbolic matrix */
347ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
3482205254eSKarl Rupp 
349a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
350a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
35125296bd5SBarry Smith 
35225296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
35325296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
35425296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
35525296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
35625296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
35725296bd5SBarry Smith   c->nonew   = 0;
3582205254eSKarl Rupp 
35925296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
36025296bd5SBarry Smith 
36125296bd5SBarry Smith   /* set MatInfo */
36225296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
36325296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
36425296bd5SBarry Smith   c->maxnz                     = ci[am];
36525296bd5SBarry Smith   c->nz                        = ci[am];
3663c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
36725296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
36825296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
36925296bd5SBarry Smith 
37025296bd5SBarry Smith #if defined(PETSC_USE_INFO)
37125296bd5SBarry Smith   if (ci[am]) {
3723c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
37325296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
37425296bd5SBarry Smith   } else {
37525296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
37625296bd5SBarry Smith   }
37725296bd5SBarry Smith #endif
37825296bd5SBarry Smith   PetscFunctionReturn(0);
37925296bd5SBarry Smith }
38025296bd5SBarry Smith 
38125296bd5SBarry Smith 
38225296bd5SBarry Smith #undef __FUNCT__
38325023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
38425023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
385e9e4536cSHong Zhang {
386e9e4536cSHong Zhang   PetscErrorCode     ierr;
387e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
388bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
38925c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
390e9e4536cSHong Zhang   MatScalar          *ca;
391e9e4536cSHong Zhang   PetscReal          afill;
392bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
3930298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
394e9e4536cSHong Zhang 
395e9e4536cSHong Zhang   PetscFunctionBegin;
396bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
397bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
398bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
399bd958071SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
400bd958071SHong Zhang   ci[0] = 0;
401bd958071SHong Zhang 
402bd958071SHong Zhang   /* create and initialize a linked list */
403bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
404bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
405bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
406bd958071SHong Zhang 
407bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
408bd958071SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
409bd958071SHong Zhang   current_space = free_space;
410bd958071SHong Zhang 
411bd958071SHong Zhang   /* Determine ci and cj */
412bd958071SHong Zhang   for (i=0; i<am; i++) {
413bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
414bd958071SHong Zhang     aj   = a->j + ai[i];
415bd958071SHong Zhang     for (j=0; j<anzi; j++) {
416bd958071SHong Zhang       brow = aj[j];
417bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
418bd958071SHong Zhang       bj   = b->j + bi[brow];
419bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
420bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
421bd958071SHong Zhang     }
422bd958071SHong Zhang     cnzi = lnk[0];
423bd958071SHong Zhang 
424bd958071SHong Zhang     /* If free space is not available, make more free space */
425bd958071SHong Zhang     /* Double the amount of total space in the list */
426bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
427bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
428bd958071SHong Zhang       ndouble++;
429bd958071SHong Zhang     }
430bd958071SHong Zhang 
431bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
432bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4332205254eSKarl Rupp 
434bd958071SHong Zhang     current_space->array           += cnzi;
435bd958071SHong Zhang     current_space->local_used      += cnzi;
436bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4372205254eSKarl Rupp 
438bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
439bd958071SHong Zhang   }
440bd958071SHong Zhang 
441bd958071SHong Zhang   /* Column indices are in the list of free space */
442bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
443bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
444bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
445bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
446bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
447e9e4536cSHong Zhang 
448e9e4536cSHong Zhang   /* Allocate space for ca */
449bd958071SHong Zhang   /*-----------------------*/
450e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
451e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
452e9e4536cSHong Zhang 
453e9e4536cSHong Zhang   /* put together the new symbolic matrix */
454ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
4552205254eSKarl Rupp 
456a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
457a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
458e9e4536cSHong Zhang 
459e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
460e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
461e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
462e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
463e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
464e9e4536cSHong Zhang   c->nonew   = 0;
4652205254eSKarl Rupp 
46625023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
467e9e4536cSHong Zhang 
468e9e4536cSHong Zhang   /* set MatInfo */
469e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
470e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
471e9e4536cSHong Zhang   c->maxnz                     = ci[am];
472e9e4536cSHong Zhang   c->nz                        = ci[am];
473bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
474e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
475e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
476e9e4536cSHong Zhang 
477e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
478e9e4536cSHong Zhang   if (ci[am]) {
479bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
480e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
481e9e4536cSHong Zhang   } else {
482e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
483e9e4536cSHong Zhang   }
484e9e4536cSHong Zhang #endif
485e9e4536cSHong Zhang   PetscFunctionReturn(0);
486e9e4536cSHong Zhang }
487e9e4536cSHong Zhang 
4880ced3a2bSJed Brown #undef __FUNCT__
4890ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
4900ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
4910ced3a2bSJed Brown {
4920ced3a2bSJed Brown   PetscErrorCode     ierr;
4930ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4940ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
4950ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
4960ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
4970ced3a2bSJed Brown   PetscReal          afill;
4980ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
4990298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
5000ced3a2bSJed Brown   PetscHeap          h;
5010ced3a2bSJed Brown 
5020ced3a2bSJed Brown   PetscFunctionBegin;
503cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
5040ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
5050ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
5060ced3a2bSJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5070ced3a2bSJed Brown   ci[0] = 0;
5080ced3a2bSJed Brown 
5090ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
5100ced3a2bSJed Brown   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5110ced3a2bSJed Brown   current_space = free_space;
5120ced3a2bSJed Brown 
5130ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
5140ced3a2bSJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
5150ced3a2bSJed Brown 
5160ced3a2bSJed Brown   /* Determine ci and cj */
5170ced3a2bSJed Brown   for (i=0; i<am; i++) {
5180ced3a2bSJed 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 */
5190ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5200ced3a2bSJed Brown     ci[i+1] = ci[i];
5210ced3a2bSJed Brown     /* Populate the min heap */
5220ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5230ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5240ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5250ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5260ced3a2bSJed Brown       }
5270ced3a2bSJed Brown     }
5280ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5290ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5300ced3a2bSJed Brown     while (j >= 0) {
5310ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
5320ced3a2bSJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
5330ced3a2bSJed Brown         ndouble++;
5340ced3a2bSJed Brown       }
5350ced3a2bSJed Brown       *(current_space->array++) = col;
5360ced3a2bSJed Brown       current_space->local_used++;
5370ced3a2bSJed Brown       current_space->local_remaining--;
5380ced3a2bSJed Brown       ci[i+1]++;
5390ced3a2bSJed Brown 
5400ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5410ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5420ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5430ced3a2bSJed Brown         PetscInt j2,col2;
5440ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5450ced3a2bSJed Brown         if (col2 != col) break;
5460ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5470ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5480ced3a2bSJed Brown       }
5490ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5500ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5510ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5520ced3a2bSJed Brown     }
5530ced3a2bSJed Brown   }
5540ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5550ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5560ced3a2bSJed Brown 
5570ced3a2bSJed Brown   /* Column indices are in the list of free space */
5580ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5590ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
5600ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5610ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5620ced3a2bSJed Brown 
5630ced3a2bSJed Brown   /* put together the new symbolic matrix */
564ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
5652205254eSKarl Rupp 
5660ced3a2bSJed Brown   (*C)->rmap->bs = A->rmap->bs;
5670ced3a2bSJed Brown   (*C)->cmap->bs = B->cmap->bs;
5680ced3a2bSJed Brown 
5690ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5700ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5710ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5720ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5730ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5740ced3a2bSJed Brown   c->nonew   = 0;
57526fbe8dcSKarl Rupp 
57689d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5770ced3a2bSJed Brown 
5780ced3a2bSJed Brown   /* set MatInfo */
5790ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5800ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5810ced3a2bSJed Brown   c->maxnz                     = ci[am];
5820ced3a2bSJed Brown   c->nz                        = ci[am];
5830ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5840ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5850ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5860ced3a2bSJed Brown 
5870ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5880ced3a2bSJed Brown   if (ci[am]) {
5890ced3a2bSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
5900ced3a2bSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
5910ced3a2bSJed Brown   } else {
5920ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
5930ced3a2bSJed Brown   }
5940ced3a2bSJed Brown #endif
5950ced3a2bSJed Brown   PetscFunctionReturn(0);
5960ced3a2bSJed Brown }
597e9e4536cSHong Zhang 
5988a07c6f1SJed Brown #undef __FUNCT__
5998a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
6008a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
6018a07c6f1SJed Brown {
6028a07c6f1SJed Brown   PetscErrorCode     ierr;
6038a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6048a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
6058a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6068a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6078a07c6f1SJed Brown   PetscReal          afill;
6088a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6090298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6108a07c6f1SJed Brown   PetscHeap          h;
6118a07c6f1SJed Brown   PetscBT            bt;
6128a07c6f1SJed Brown 
6138a07c6f1SJed Brown   PetscFunctionBegin;
614cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6158a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6168a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
6178a07c6f1SJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6188a07c6f1SJed Brown   ci[0] = 0;
6198a07c6f1SJed Brown 
6208a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
6218a07c6f1SJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6222205254eSKarl Rupp 
6238a07c6f1SJed Brown   current_space = free_space;
6248a07c6f1SJed Brown 
6258a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
6268a07c6f1SJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
6278a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6288a07c6f1SJed Brown 
6298a07c6f1SJed Brown   /* Determine ci and cj */
6308a07c6f1SJed Brown   for (i=0; i<am; i++) {
6318a07c6f1SJed 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 */
6328a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6338a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6348a07c6f1SJed Brown     ci[i+1] = ci[i];
6358a07c6f1SJed Brown     /* Populate the min heap */
6368a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6378a07c6f1SJed Brown       PetscInt brow = acol[j];
6388a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6398a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6408a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6418a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6428a07c6f1SJed Brown           bb[j]++;
6438a07c6f1SJed Brown           break;
6448a07c6f1SJed Brown         }
6458a07c6f1SJed Brown       }
6468a07c6f1SJed Brown     }
6478a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6488a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6498a07c6f1SJed Brown     while (j >= 0) {
6508a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6510298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
6528a07c6f1SJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
6538a07c6f1SJed Brown         ndouble++;
6548a07c6f1SJed Brown       }
6558a07c6f1SJed Brown       *(current_space->array++) = col;
6568a07c6f1SJed Brown       current_space->local_used++;
6578a07c6f1SJed Brown       current_space->local_remaining--;
6588a07c6f1SJed Brown       ci[i+1]++;
6598a07c6f1SJed Brown 
6608a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6618a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6628a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6638a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6648a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6658a07c6f1SJed Brown           bb[j]++;
6668a07c6f1SJed Brown           break;
6678a07c6f1SJed Brown         }
6688a07c6f1SJed Brown       }
6698a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6708a07c6f1SJed Brown     }
6718a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6728a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6738a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6748a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6758a07c6f1SJed Brown     }
6768a07c6f1SJed Brown   }
6778a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6788a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6798a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6808a07c6f1SJed Brown 
6818a07c6f1SJed Brown   /* Column indices are in the list of free space */
6828a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6838a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
6848a07c6f1SJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6858a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6868a07c6f1SJed Brown 
6878a07c6f1SJed Brown   /* put together the new symbolic matrix */
688ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
6892205254eSKarl Rupp 
6908a07c6f1SJed Brown   (*C)->rmap->bs = A->rmap->bs;
6918a07c6f1SJed Brown   (*C)->cmap->bs = B->cmap->bs;
6928a07c6f1SJed Brown 
6938a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6948a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6958a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
6968a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
6978a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
6988a07c6f1SJed Brown   c->nonew   = 0;
69926fbe8dcSKarl Rupp 
70089d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
7018a07c6f1SJed Brown 
7028a07c6f1SJed Brown   /* set MatInfo */
7038a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
7048a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
7058a07c6f1SJed Brown   c->maxnz                     = ci[am];
7068a07c6f1SJed Brown   c->nz                        = ci[am];
7078a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7088a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7098a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7108a07c6f1SJed Brown 
7118a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7128a07c6f1SJed Brown   if (ci[am]) {
7138a07c6f1SJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
7148a07c6f1SJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
7158a07c6f1SJed Brown   } else {
7168a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7178a07c6f1SJed Brown   }
7188a07c6f1SJed Brown #endif
7198a07c6f1SJed Brown   PetscFunctionReturn(0);
7208a07c6f1SJed Brown }
7218a07c6f1SJed Brown 
722cd093f1dSJed Brown #undef __FUNCT__
72358cf0668SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
724cd093f1dSJed Brown /* concatenate unique entries and then sort */
72558cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
726cd093f1dSJed Brown {
727cd093f1dSJed Brown   PetscErrorCode     ierr;
728cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
729cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
730cd093f1dSJed Brown   PetscInt           *ci,*cj;
731cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
732cd093f1dSJed Brown   PetscReal          afill;
733cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
734cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
735cd093f1dSJed Brown   char               *seen;
736cd093f1dSJed Brown 
737cd093f1dSJed Brown   PetscFunctionBegin;
738cd093f1dSJed Brown   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
739cd093f1dSJed Brown   ci[0] = 0;
740cd093f1dSJed Brown 
741cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
742cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
743cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
744cd093f1dSJed Brown   ierr = PetscMalloc(bn*sizeof(char),&seen);CHKERRQ(ierr);
745cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
746cd093f1dSJed Brown 
747cd093f1dSJed Brown   /* Determine ci and cj */
748cd093f1dSJed Brown   for (i=0; i<am; i++) {
749cd093f1dSJed Brown     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
750cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
751cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
752cd093f1dSJed Brown     /* Pack segrow */
753cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
754cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
755cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
756cd093f1dSJed Brown         PetscInt bcol = bj[k];
757cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
758cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
759cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
760cd093f1dSJed Brown           *slot = bcol;
761cd093f1dSJed Brown           seen[bcol] = 1;
762cd093f1dSJed Brown           packlen++;
763cd093f1dSJed Brown         }
764cd093f1dSJed Brown       }
765cd093f1dSJed Brown     }
766cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
767cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
768cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
769cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
770cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
771cd093f1dSJed Brown   }
772cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
773cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
774cd093f1dSJed Brown 
775cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
776cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
777cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
778cd093f1dSJed Brown 
779cd093f1dSJed Brown   /* put together the new symbolic matrix */
780cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
781cd093f1dSJed Brown 
782cd093f1dSJed Brown   (*C)->rmap->bs = A->rmap->bs;
783cd093f1dSJed Brown   (*C)->cmap->bs = B->cmap->bs;
784cd093f1dSJed Brown 
785cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
786cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
787cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
788cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
789cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
790cd093f1dSJed Brown   c->nonew   = 0;
791cd093f1dSJed Brown 
792cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
793cd093f1dSJed Brown 
794cd093f1dSJed Brown   /* set MatInfo */
795cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
796cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
797cd093f1dSJed Brown   c->maxnz                     = ci[am];
798cd093f1dSJed Brown   c->nz                        = ci[am];
799cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
800cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
801cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
802cd093f1dSJed Brown 
803cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
804cd093f1dSJed Brown   if (ci[am]) {
805cd093f1dSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
806cd093f1dSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
807cd093f1dSJed Brown   } else {
808cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
809cd093f1dSJed Brown   }
810cd093f1dSJed Brown #endif
811cd093f1dSJed Brown   PetscFunctionReturn(0);
812cd093f1dSJed Brown }
813cd093f1dSJed Brown 
814d2853540SHong Zhang /* This routine is not used. Should be removed! */
815bc011b1eSHong Zhang #undef __FUNCT__
8166fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
8176fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8185df89d91SHong Zhang {
819bc011b1eSHong Zhang   PetscErrorCode ierr;
820bc011b1eSHong Zhang 
821bc011b1eSHong Zhang   PetscFunctionBegin;
822bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
8233ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
8246fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8253ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
826bc011b1eSHong Zhang   }
8273ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
8286fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8293ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
830bc011b1eSHong Zhang   PetscFunctionReturn(0);
831bc011b1eSHong Zhang }
832bc011b1eSHong Zhang 
833bc011b1eSHong Zhang #undef __FUNCT__
8342128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
8352128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8362128a86cSHong Zhang {
8372128a86cSHong Zhang   PetscErrorCode      ierr;
8384c7df5ccSHong Zhang   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
83940192850SHong Zhang   Mat_MatMatTransMult *abt=a->abt;
8402128a86cSHong Zhang 
8412128a86cSHong Zhang   PetscFunctionBegin;
84240192850SHong Zhang   ierr = (abt->destroy)(A);CHKERRQ(ierr);
84340192850SHong Zhang   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
84440192850SHong Zhang   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
84540192850SHong Zhang   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
84640192850SHong Zhang   ierr = PetscFree(abt);CHKERRQ(ierr);
8472128a86cSHong Zhang   PetscFunctionReturn(0);
8482128a86cSHong Zhang }
8492128a86cSHong Zhang 
8502128a86cSHong Zhang #undef __FUNCT__
8516fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
8526fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
853bc011b1eSHong Zhang {
8545df89d91SHong Zhang   PetscErrorCode      ierr;
85581d82fe4SHong Zhang   Mat                 Bt;
85681d82fe4SHong Zhang   PetscInt            *bti,*btj;
85740192850SHong Zhang   Mat_MatMatTransMult *abt;
8584c7df5ccSHong Zhang   Mat_SeqAIJ          *c;
859d2853540SHong Zhang 
86081d82fe4SHong Zhang   PetscFunctionBegin;
86181d82fe4SHong Zhang   /* create symbolic Bt */
86281d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8630298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
8642205254eSKarl Rupp 
865c3e5699bSSatish Balay   Bt->rmap->bs = A->cmap->bs;
866c3e5699bSSatish Balay   Bt->cmap->bs = B->cmap->bs;
86781d82fe4SHong Zhang 
86881d82fe4SHong Zhang   /* get symbolic C=A*Bt */
86981d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
87081d82fe4SHong Zhang 
8712128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
87240192850SHong Zhang   ierr   = PetscNew(Mat_MatMatTransMult,&abt);CHKERRQ(ierr);
8734c7df5ccSHong Zhang   c      = (Mat_SeqAIJ*)(*C)->data;
87440192850SHong Zhang   c->abt = abt;
8752128a86cSHong Zhang 
87640192850SHong Zhang   abt->usecoloring = PETSC_FALSE;
87740192850SHong Zhang   abt->destroy     = (*C)->ops->destroy;
8782128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8792128a86cSHong Zhang 
88040192850SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr);
88140192850SHong Zhang   if (abt->usecoloring) {
882b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
883b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
884*335efc43SPeter Brune     MatColoring          coloring;
8852128a86cSHong Zhang     ISColoring           iscoloring;
8862128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
8874d478ae7SHong Zhang     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
8884d478ae7SHong Zhang     /* inode causes memory problem, don't know why */
8894d478ae7SHong Zhang     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
890e8354b3eSHong Zhang 
891*335efc43SPeter Brune     ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr);
892*335efc43SPeter Brune     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
893*335efc43SPeter Brune     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
894*335efc43SPeter Brune     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
895*335efc43SPeter Brune     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
896*335efc43SPeter Brune     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
897b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
8982205254eSKarl Rupp 
89940192850SHong Zhang     abt->matcoloring = matcoloring;
9002205254eSKarl Rupp 
9012128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
9022128a86cSHong Zhang 
9032128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
9042128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
9052128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9062128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
9070298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9082205254eSKarl Rupp 
9092128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
91040192850SHong Zhang     abt->Bt_den   = Bt_dense;
9112128a86cSHong Zhang 
9122128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9132128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9142128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9150298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9162205254eSKarl Rupp 
9172128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
91840192850SHong Zhang     abt->ABt_den  = C_dense;
919f94ccd6cSHong Zhang 
920f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9211ce71dffSSatish Balay     {
922f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
923c40ebe3bSHong Zhang       ierr = PetscInfo7(*C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr);
9241ce71dffSSatish Balay     }
925f94ccd6cSHong Zhang #endif
9262128a86cSHong Zhang   }
927e99be685SHong Zhang   /* clean up */
928e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
929e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9305df89d91SHong Zhang   PetscFunctionReturn(0);
9315df89d91SHong Zhang }
9325df89d91SHong Zhang 
9335df89d91SHong Zhang #undef __FUNCT__
9346fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9356fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9365df89d91SHong Zhang {
9375df89d91SHong Zhang   PetscErrorCode      ierr;
9385df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
939e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
94081d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9415df89d91SHong Zhang   PetscLogDouble      flops=0.0;
942aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
94340192850SHong Zhang   Mat_MatMatTransMult *abt = c->abt;
9445df89d91SHong Zhang 
9455df89d91SHong Zhang   PetscFunctionBegin;
94658ed3ceaSHong Zhang   /* clear old values in C */
94758ed3ceaSHong Zhang   if (!c->a) {
94858ed3ceaSHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
94958ed3ceaSHong Zhang     c->a      = ca;
95058ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
95158ed3ceaSHong Zhang   } else {
95258ed3ceaSHong Zhang     ca =  c->a;
95358ed3ceaSHong Zhang   }
95458ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
95558ed3ceaSHong Zhang 
95640192850SHong Zhang   if (abt->usecoloring) {
95740192850SHong Zhang     MatTransposeColoring matcoloring = abt->matcoloring;
95840192850SHong Zhang     Mat                  Bt_dense,C_dense = abt->ABt_den;
959c8db22f6SHong Zhang 
960b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
96140192850SHong Zhang     Bt_dense = abt->Bt_den;
962b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
963c8db22f6SHong Zhang 
964c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9652128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
966c8db22f6SHong Zhang 
9672128a86cSHong Zhang     /* Recover C from C_dense */
968b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
969c8db22f6SHong Zhang     PetscFunctionReturn(0);
970c8db22f6SHong Zhang   }
971c8db22f6SHong Zhang 
9721fa1209cSHong Zhang   for (i=0; i<cm; i++) {
97381d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
97481d82fe4SHong Zhang     acol = aj + ai[i];
9756973769bSHong Zhang     aval = aa + ai[i];
9761fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9771fa1209cSHong Zhang     ccol = cj + ci[i];
9786973769bSHong Zhang     cval = ca + ci[i];
9791fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
98081d82fe4SHong Zhang       brow = ccol[j];
98181d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
98281d82fe4SHong Zhang       bcol = bj + bi[brow];
9836973769bSHong Zhang       bval = ba + bi[brow];
9846973769bSHong Zhang 
98581d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
98681d82fe4SHong Zhang       nexta = 0; nextb = 0;
98781d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
9887b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
98981d82fe4SHong Zhang         if (nexta == anzi) break;
9907b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
99181d82fe4SHong Zhang         if (nextb == bnzj) break;
99281d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
9936973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
99481d82fe4SHong Zhang           nexta++; nextb++;
99581d82fe4SHong Zhang           flops += 2;
9961fa1209cSHong Zhang         }
9971fa1209cSHong Zhang       }
99881d82fe4SHong Zhang     }
99981d82fe4SHong Zhang   }
100081d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100181d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100281d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10035df89d91SHong Zhang   PetscFunctionReturn(0);
10045df89d91SHong Zhang }
10055df89d91SHong Zhang 
10065df89d91SHong Zhang #undef __FUNCT__
100775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
10080adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10090adebc6cSBarry Smith {
10105df89d91SHong Zhang   PetscErrorCode ierr;
10115df89d91SHong Zhang 
10125df89d91SHong Zhang   PetscFunctionBegin;
10135df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
101407706670SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
101575648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
101607706670SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
10175df89d91SHong Zhang   }
101807706670SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
101975648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
102007706670SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
10215df89d91SHong Zhang   PetscFunctionReturn(0);
10225df89d91SHong Zhang }
10235df89d91SHong Zhang 
10245df89d91SHong Zhang #undef __FUNCT__
102575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
102675648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10275df89d91SHong Zhang {
1028bc011b1eSHong Zhang   PetscErrorCode ierr;
1029bc011b1eSHong Zhang   Mat            At;
103038baddfdSBarry Smith   PetscInt       *ati,*atj;
1031bc011b1eSHong Zhang 
1032bc011b1eSHong Zhang   PetscFunctionBegin;
1033bc011b1eSHong Zhang   /* create symbolic At */
1034bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10350298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
10362205254eSKarl Rupp 
1037a2f3521dSMark F. Adams   At->rmap->bs = A->cmap->bs;
1038a2f3521dSMark F. Adams   At->cmap->bs = B->cmap->bs;
1039bc011b1eSHong Zhang 
1040bc011b1eSHong Zhang   /* get symbolic C=At*B */
1041bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1042bc011b1eSHong Zhang 
1043bc011b1eSHong Zhang   /* clean up */
10446bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1045bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1046bc011b1eSHong Zhang   PetscFunctionReturn(0);
1047bc011b1eSHong Zhang }
1048bc011b1eSHong Zhang 
1049bc011b1eSHong Zhang #undef __FUNCT__
105075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
105175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1052bc011b1eSHong Zhang {
1053bc011b1eSHong Zhang   PetscErrorCode ierr;
10540fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1055d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1056d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1057d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1058aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1059bc011b1eSHong Zhang 
1060bc011b1eSHong Zhang   PetscFunctionBegin;
1061aa1aec99SHong Zhang   if (!c->a) {
1062aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
10632205254eSKarl Rupp 
1064aa1aec99SHong Zhang     c->a      = ca;
1065aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1066aa1aec99SHong Zhang   } else {
1067aa1aec99SHong Zhang     ca = c->a;
1068aa1aec99SHong Zhang   }
1069bc011b1eSHong Zhang   /* clear old values in C */
1070bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1071bc011b1eSHong Zhang 
1072bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1073bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1074bc011b1eSHong Zhang     bj   = b->j + bi[i];
1075bc011b1eSHong Zhang     ba   = b->a + bi[i];
1076bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1077bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1078bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1079bc011b1eSHong Zhang       nextb = 0;
10800fbc74f4SHong Zhang       crow  = *aj++;
1081bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1082bc011b1eSHong Zhang       caj   = ca + ci[crow];
1083bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1084bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
10850fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
10860fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1087bc011b1eSHong Zhang           nextb++;
1088bc011b1eSHong Zhang         }
1089bc011b1eSHong Zhang       }
1090bc011b1eSHong Zhang       flops += 2*bnzi;
10910fbc74f4SHong Zhang       aa++;
1092bc011b1eSHong Zhang     }
1093bc011b1eSHong Zhang   }
1094bc011b1eSHong Zhang 
1095bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1096bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1097bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1098bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1099bc011b1eSHong Zhang   PetscFunctionReturn(0);
1100bc011b1eSHong Zhang }
11019123193aSHong Zhang 
11029123193aSHong Zhang #undef __FUNCT__
11039123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
11049123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11059123193aSHong Zhang {
11069123193aSHong Zhang   PetscErrorCode ierr;
11079123193aSHong Zhang 
11089123193aSHong Zhang   PetscFunctionBegin;
11099123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
11103ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11119123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11123ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11139123193aSHong Zhang   }
11143ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11159123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11164614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11179123193aSHong Zhang   PetscFunctionReturn(0);
11189123193aSHong Zhang }
1119edd81eecSMatthew Knepley 
11209123193aSHong Zhang #undef __FUNCT__
11219123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11229123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11239123193aSHong Zhang {
11249123193aSHong Zhang   PetscErrorCode ierr;
11259123193aSHong Zhang 
11269123193aSHong Zhang   PetscFunctionBegin;
11275a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11282205254eSKarl Rupp 
1129d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11309123193aSHong Zhang   PetscFunctionReturn(0);
11319123193aSHong Zhang }
11329123193aSHong Zhang 
11339123193aSHong Zhang #undef __FUNCT__
11349123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11359123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11369123193aSHong Zhang {
1137f32d5d43SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
11389123193aSHong Zhang   PetscErrorCode    ierr;
1139daf57211SHong Zhang   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1140daf57211SHong Zhang   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1141daf57211SHong Zhang   const PetscInt    *aj;
1142daf57211SHong Zhang   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=B->rmap->n,am=A->rmap->n;
1143daf57211SHong Zhang   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
11449123193aSHong Zhang 
11459123193aSHong Zhang   PetscFunctionBegin;
1146f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1147e32f2f54SBarry 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);
1148e32f2f54SBarry 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);
1149e32f2f54SBarry 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);
11508c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
11518c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1152f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1153730858b9SHong Zhang   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1154f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1155f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1156f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1157f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1158f32d5d43SBarry Smith       aj = a->j + a->i[i];
1159f32d5d43SBarry Smith       aa = a->a + a->i[i];
1160f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1161730858b9SHong Zhang         aatmp = aa[j]; ajtmp = aj[j];
1162730858b9SHong Zhang         r1 += aatmp*b1[ajtmp];
1163730858b9SHong Zhang         r2 += aatmp*b2[ajtmp];
1164730858b9SHong Zhang         r3 += aatmp*b3[ajtmp];
1165730858b9SHong Zhang         r4 += aatmp*b4[ajtmp];
11669123193aSHong Zhang       }
1167730858b9SHong Zhang       c1[i] = r1;
1168730858b9SHong Zhang       c2[i] = r2;
1169730858b9SHong Zhang       c3[i] = r3;
1170730858b9SHong Zhang       c4[i] = r4;
1171f32d5d43SBarry Smith     }
1172730858b9SHong Zhang     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1173730858b9SHong Zhang     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1174f32d5d43SBarry Smith   }
1175f32d5d43SBarry Smith   for (; col<cn; col++) {   /* over extra columns of C */
1176f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1177f32d5d43SBarry Smith       r1 = 0.0;
1178f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1179f32d5d43SBarry Smith       aj = a->j + a->i[i];
1180f32d5d43SBarry Smith       aa = a->a + a->i[i];
1181f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1182730858b9SHong Zhang         r1 += aa[j]*b1[aj[j]];
1183f32d5d43SBarry Smith       }
1184730858b9SHong Zhang       c1[i] = r1;
1185f32d5d43SBarry Smith     }
1186f32d5d43SBarry Smith     b1 += bm;
1187730858b9SHong Zhang     c1 += am;
1188f32d5d43SBarry Smith   }
1189dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
11908c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
11918c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1192f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1193f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1194f32d5d43SBarry Smith   PetscFunctionReturn(0);
1195f32d5d43SBarry Smith }
1196f32d5d43SBarry Smith 
1197f32d5d43SBarry Smith /*
11984324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1199f32d5d43SBarry Smith */
1200f32d5d43SBarry Smith #undef __FUNCT__
1201f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1202f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1203f32d5d43SBarry Smith {
1204f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1205f32d5d43SBarry Smith   PetscErrorCode ierr;
1206dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1207dd6ea824SBarry Smith   MatScalar      *aa;
1208d0f46423SBarry 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;
12094324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1210f32d5d43SBarry Smith 
1211f32d5d43SBarry Smith   PetscFunctionBegin;
1212f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
12138c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12148c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1215f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12164324174eSBarry Smith 
12174324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12184324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12194324174eSBarry Smith       colam = col*am;
12204324174eSBarry Smith       arm   = a->compressedrow.nrows;
12214324174eSBarry Smith       ii    = a->compressedrow.i;
12224324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12234324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12244324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12254324174eSBarry Smith         n  = ii[i+1] - ii[i];
12264324174eSBarry Smith         aj = a->j + ii[i];
12274324174eSBarry Smith         aa = a->a + ii[i];
12284324174eSBarry Smith         for (j=0; j<n; j++) {
12294324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12304324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12314324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12324324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12334324174eSBarry Smith         }
12344324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12354324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12364324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12374324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12384324174eSBarry Smith       }
12394324174eSBarry Smith       b1 += bm4;
12404324174eSBarry Smith       b2 += bm4;
12414324174eSBarry Smith       b3 += bm4;
12424324174eSBarry Smith       b4 += bm4;
12434324174eSBarry Smith     }
12444324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12454324174eSBarry Smith       colam = col*am;
12464324174eSBarry Smith       arm   = a->compressedrow.nrows;
12474324174eSBarry Smith       ii    = a->compressedrow.i;
12484324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12494324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12504324174eSBarry Smith         r1 = 0.0;
12514324174eSBarry Smith         n  = ii[i+1] - ii[i];
12524324174eSBarry Smith         aj = a->j + ii[i];
12534324174eSBarry Smith         aa = a->a + ii[i];
12544324174eSBarry Smith 
12554324174eSBarry Smith         for (j=0; j<n; j++) {
12564324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
12574324174eSBarry Smith         }
1258a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
12594324174eSBarry Smith       }
12604324174eSBarry Smith       b1 += bm;
12614324174eSBarry Smith     }
12624324174eSBarry Smith   } else {
1263f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1264f32d5d43SBarry Smith       colam = col*am;
1265f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1266f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1267f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1268f32d5d43SBarry Smith         aj = a->j + a->i[i];
1269f32d5d43SBarry Smith         aa = a->a + a->i[i];
1270f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1271f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1272f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1273f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1274f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1275f32d5d43SBarry Smith         }
1276f32d5d43SBarry Smith         c[colam + i]       += r1;
1277f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1278f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1279f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1280f32d5d43SBarry Smith       }
1281f32d5d43SBarry Smith       b1 += bm4;
1282f32d5d43SBarry Smith       b2 += bm4;
1283f32d5d43SBarry Smith       b3 += bm4;
1284f32d5d43SBarry Smith       b4 += bm4;
1285f32d5d43SBarry Smith     }
1286f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1287a2ea699eSBarry Smith       colam = col*am;
1288f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1289f32d5d43SBarry Smith         r1 = 0.0;
1290f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1291f32d5d43SBarry Smith         aj = a->j + a->i[i];
1292f32d5d43SBarry Smith         aa = a->a + a->i[i];
1293f32d5d43SBarry Smith 
1294f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1295f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1296f32d5d43SBarry Smith         }
1297a2ea699eSBarry Smith         c[colam + i] += r1;
1298f32d5d43SBarry Smith       }
1299f32d5d43SBarry Smith       b1 += bm;
1300f32d5d43SBarry Smith     }
13014324174eSBarry Smith   }
1302dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13038c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
13048c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13059123193aSHong Zhang   PetscFunctionReturn(0);
13069123193aSHong Zhang }
1307b1683b59SHong Zhang 
1308b1683b59SHong Zhang #undef __FUNCT__
1309b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1310b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1311c8db22f6SHong Zhang {
1312c8db22f6SHong Zhang   PetscErrorCode ierr;
13132f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13142f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13152f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13162f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13172f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13182f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1319c8db22f6SHong Zhang 
1320c8db22f6SHong Zhang   PetscFunctionBegin;
13212f5f1f90SHong Zhang   btval_den=btdense->v;
13222f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13232f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13242f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13252f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1326d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13272f5f1f90SHong Zhang       btcol = bj + bi[col];
13282f5f1f90SHong Zhang       btval = ba + bi[col];
13292f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1330d2853540SHong Zhang       for (j=0; j<anz; j++) {
13312f5f1f90SHong Zhang         brow            = btcol[j];
13322f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1333c8db22f6SHong Zhang       }
1334c8db22f6SHong Zhang     }
13352f5f1f90SHong Zhang     btval_den += m;
1336c8db22f6SHong Zhang   }
1337c8db22f6SHong Zhang   PetscFunctionReturn(0);
1338c8db22f6SHong Zhang }
1339c8db22f6SHong Zhang 
1340c8db22f6SHong Zhang #undef __FUNCT__
1341b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1342b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13438972f759SHong Zhang {
13448972f759SHong Zhang   PetscErrorCode ierr;
1345b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
1346077f23c2SHong Zhang   PetscScalar    *ca_den,*ca_den_ptr,*ca=csp->a;
1347f99a636bSHong Zhang   PetscInt       k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1348e88777f2SHong Zhang   PetscInt       brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1349077f23c2SHong Zhang   PetscInt       nrows,*row,*idx;
1350077f23c2SHong Zhang   PetscInt       *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
13518972f759SHong Zhang 
13528972f759SHong Zhang   PetscFunctionBegin;
1353a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1354f99a636bSHong Zhang 
1355077f23c2SHong Zhang   if (brows > 0) {
1356077f23c2SHong Zhang     PetscInt *lstart,row_end,row_start;
1357980ae229SHong Zhang     lstart = matcoloring->lstart;
1358eeb4fd02SHong Zhang     ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr);
1359eeb4fd02SHong Zhang 
1360077f23c2SHong Zhang     row_end = brows;
1361eeb4fd02SHong Zhang     if (row_end > m) row_end = m;
1362077f23c2SHong Zhang     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1363077f23c2SHong Zhang       ca_den_ptr = ca_den;
1364980ae229SHong Zhang       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1365eeb4fd02SHong Zhang         nrows = matcoloring->nrows[k];
1366eeb4fd02SHong Zhang         row   = rows  + colorforrow[k];
1367eeb4fd02SHong Zhang         idx   = den2sp + colorforrow[k];
1368eeb4fd02SHong Zhang         for (l=lstart[k]; l<nrows; l++) {
1369eeb4fd02SHong Zhang           if (row[l] >= row_end) {
1370eeb4fd02SHong Zhang             lstart[k] = l;
1371eeb4fd02SHong Zhang             break;
1372eeb4fd02SHong Zhang           } else {
1373077f23c2SHong Zhang             ca[idx[l]] = ca_den_ptr[row[l]];
1374eeb4fd02SHong Zhang           }
1375eeb4fd02SHong Zhang         }
1376077f23c2SHong Zhang         ca_den_ptr += m;
1377eeb4fd02SHong Zhang       }
1378077f23c2SHong Zhang       row_end += brows;
1379eeb4fd02SHong Zhang       if (row_end > m) row_end = m;
1380eeb4fd02SHong Zhang     }
1381077f23c2SHong Zhang   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1382077f23c2SHong Zhang     ca_den_ptr = ca_den;
1383077f23c2SHong Zhang     for (k=0; k<ncolors; k++) {
1384077f23c2SHong Zhang       nrows = matcoloring->nrows[k];
1385077f23c2SHong Zhang       row   = rows  + colorforrow[k];
1386077f23c2SHong Zhang       idx   = den2sp + colorforrow[k];
1387077f23c2SHong Zhang       for (l=0; l<nrows; l++) {
1388077f23c2SHong Zhang         ca[idx[l]] = ca_den_ptr[row[l]];
1389077f23c2SHong Zhang       }
1390077f23c2SHong Zhang       ca_den_ptr += m;
1391077f23c2SHong Zhang     }
1392f99a636bSHong Zhang   }
1393f99a636bSHong Zhang 
1394a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
1395f99a636bSHong Zhang #if defined(PETSC_USE_INFO)
1396077f23c2SHong Zhang   if (matcoloring->brows > 0) {
1397f99a636bSHong Zhang     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1398e88777f2SHong Zhang   } else {
1399077f23c2SHong Zhang     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1400e88777f2SHong Zhang   }
1401f99a636bSHong Zhang #endif
14028972f759SHong Zhang   PetscFunctionReturn(0);
14038972f759SHong Zhang }
14048972f759SHong Zhang 
1405e99be685SHong Zhang /*
1406e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1407e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1408e88777f2SHong Zhang  spidx[], index of a->a, to be used in MatTransposeColoringCreate_SeqAIJ().
1409e99be685SHong Zhang  */
1410e99be685SHong Zhang #undef __FUNCT__
1411e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
14121a83f524SJed 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)
1413e99be685SHong Zhang {
1414e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1415e99be685SHong Zhang   PetscErrorCode ierr;
1416e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1417e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1418e99be685SHong Zhang   PetscInt       *cspidx;
1419e99be685SHong Zhang 
1420e99be685SHong Zhang   PetscFunctionBegin;
1421e99be685SHong Zhang   *nn = n;
1422e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1423e99be685SHong Zhang   if (symmetric) {
1424ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
14251a83f524SJed Brown     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr);
1426e99be685SHong Zhang   } else {
1427e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1428e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1429e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1430e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1431e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1432e99be685SHong Zhang     jj   = a->j;
1433e99be685SHong Zhang     for (i=0; i<nz; i++) {
1434e99be685SHong Zhang       collengths[jj[i]]++;
1435e99be685SHong Zhang     }
1436e99be685SHong Zhang     cia[0] = oshift;
1437e99be685SHong Zhang     for (i=0; i<n; i++) {
1438e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1439e99be685SHong Zhang     }
1440e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1441e99be685SHong Zhang     jj   = a->j;
1442e99be685SHong Zhang     for (row=0; row<m; row++) {
1443e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1444e99be685SHong Zhang       for (i=0; i<mr; i++) {
1445e99be685SHong Zhang         col = *jj++;
14462205254eSKarl Rupp 
1447e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1448e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift]  = row + oshift;
1449e99be685SHong Zhang       }
1450e99be685SHong Zhang     }
1451e99be685SHong Zhang     ierr   = PetscFree(collengths);CHKERRQ(ierr);
1452e99be685SHong Zhang     *ia    = cia; *ja = cja;
1453e99be685SHong Zhang     *spidx = cspidx;
1454e99be685SHong Zhang   }
1455e99be685SHong Zhang   PetscFunctionReturn(0);
1456e99be685SHong Zhang }
1457e99be685SHong Zhang 
1458e99be685SHong Zhang #undef __FUNCT__
1459e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
14601a83f524SJed 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)
1461e99be685SHong Zhang {
1462e99be685SHong Zhang   PetscErrorCode ierr;
1463e99be685SHong Zhang 
1464e99be685SHong Zhang   PetscFunctionBegin;
1465e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1466e99be685SHong Zhang 
1467e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1468e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1469e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1470e99be685SHong Zhang   PetscFunctionReturn(0);
1471e99be685SHong Zhang }
1472e99be685SHong Zhang 
14738972f759SHong Zhang #undef __FUNCT__
1474b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1475b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1476b1683b59SHong Zhang {
1477b1683b59SHong Zhang   PetscErrorCode ierr;
1478e88777f2SHong Zhang   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
14791a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1480b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1481b1683b59SHong Zhang   IS             *isa;
1482b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1483e88777f2SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1484e88777f2SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i,brows;
1485e88777f2SHong Zhang   PetscBool      flg;
1486b1683b59SHong Zhang 
1487b1683b59SHong Zhang   PetscFunctionBegin;
1488b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1489e99be685SHong Zhang 
14907ecccc15SHong Zhang   /* bs >1 is not being tested yet! */
1491e88777f2SHong Zhang   Nbs       = mat->cmap->N/bs;
1492b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1493e88777f2SHong Zhang   c->N      = Nbs;
1494e88777f2SHong Zhang   c->m      = c->M;
1495b1683b59SHong Zhang   c->rstart = 0;
1496e88777f2SHong Zhang   c->brows  = 100;
1497b1683b59SHong Zhang 
1498b1683b59SHong Zhang   c->ncolors = nis;
1499077f23c2SHong Zhang   ierr = PetscMalloc3(nis,PetscInt,&c->ncolumns,nis,PetscInt,&c->nrows,nis+1,PetscInt,&colorforrow);CHKERRQ(ierr);
1500e88777f2SHong Zhang   ierr = PetscMalloc((csp->nz+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr);
1501077f23c2SHong Zhang   ierr = PetscMalloc((csp->nz+1)*sizeof(PetscInt),&den2sp);CHKERRQ(ierr);
1502e88777f2SHong Zhang 
1503e88777f2SHong Zhang   brows = c->brows;
1504e88777f2SHong Zhang   ierr = PetscOptionsGetInt(NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1505e88777f2SHong Zhang   if (flg) c->brows = brows;
1506eeb4fd02SHong Zhang   if (brows > 0) {
1507980ae229SHong Zhang     ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&c->lstart);CHKERRQ(ierr);
1508e88777f2SHong Zhang   }
15092205254eSKarl Rupp 
1510d2853540SHong Zhang   colorforrow[0] = 0;
1511d2853540SHong Zhang   rows_i         = rows;
1512f99a636bSHong Zhang   den2sp_i       = den2sp;
1513b1683b59SHong Zhang 
1514d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1515e88777f2SHong Zhang   ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
15162205254eSKarl Rupp 
1517d2853540SHong Zhang   colorforcol[0] = 0;
1518d2853540SHong Zhang   columns_i      = columns;
1519d2853540SHong Zhang 
1520077f23c2SHong Zhang   /* get column-wise storage of mat */
1521077f23c2SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1522b1683b59SHong Zhang 
1523b28d80bdSHong Zhang   cm   = c->m;
1524b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1525e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1526f99a636bSHong Zhang   for (i=0; i<nis; i++) { /* loop over color */
1527b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1528b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
15292205254eSKarl Rupp 
1530b1683b59SHong Zhang     c->ncolumns[i] = n;
1531b1683b59SHong Zhang     if (n) {
1532d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1533b1683b59SHong Zhang     }
1534d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1535d2853540SHong Zhang     columns_i       += n;
1536b1683b59SHong Zhang 
1537b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1538e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1539e99be685SHong Zhang 
1540b1683b59SHong Zhang     /* loop over columns*/
1541b1683b59SHong Zhang     for (j=0; j<n; j++) {
1542b1683b59SHong Zhang       col     = is[j];
1543d2853540SHong Zhang       row_idx = cj + ci[col];
1544b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1545b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1546b1683b59SHong Zhang       for (k=0; k<m; k++) {
1547e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1548d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1549b1683b59SHong Zhang       }
1550b1683b59SHong Zhang     }
1551b1683b59SHong Zhang     /* count the number of hits */
1552b1683b59SHong Zhang     nrows = 0;
1553e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1554b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1555b1683b59SHong Zhang     }
1556b1683b59SHong Zhang     c->nrows[i]      = nrows;
1557d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1558d2853540SHong Zhang 
1559b1683b59SHong Zhang     nrows = 0;
1560e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1561b1683b59SHong Zhang       if (rowhit[j]) {
1562d2853540SHong Zhang         rows_i[nrows]   = j;
156312b89a43SHong Zhang         den2sp_i[nrows] = idxhit[j];
1564b1683b59SHong Zhang         nrows++;
1565b1683b59SHong Zhang       }
1566b1683b59SHong Zhang     }
1567e88777f2SHong Zhang     den2sp_i += nrows;
1568077f23c2SHong Zhang 
1569b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1570f99a636bSHong Zhang     rows_i += nrows;
1571b1683b59SHong Zhang   }
15720298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1573b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1574b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1575d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1576b28d80bdSHong Zhang 
1577d2853540SHong Zhang   c->colorforrow = colorforrow;
1578d2853540SHong Zhang   c->rows        = rows;
1579f99a636bSHong Zhang   c->den2sp      = den2sp;
1580d2853540SHong Zhang   c->colorforcol = colorforcol;
1581d2853540SHong Zhang   c->columns     = columns;
15822205254eSKarl Rupp 
1583f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1584b1683b59SHong Zhang   PetscFunctionReturn(0);
1585b1683b59SHong Zhang }
1586