xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 3ff4c91cb900f8ad58d90ce5eedf1ebe891d7c8c)
1be1d678aSKris Buschelman 
2d50806bdSBarry Smith /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4d50806bdSBarry Smith           C = A * B
5d50806bdSBarry Smith */
6d50806bdSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
90ced3a2bSJed Brown #include <../src/mat/utils/petscheap.h>
10c6db04a5SJed Brown #include <petscbt.h>
1107475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h>
127bab7c10SHong Zhang 
136284ec50SHong Zhang #undef __FUNCT__
145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1638baddfdSBarry Smith {
17dfbe8321SBarry Smith   PetscErrorCode ierr;
188a07c6f1SJed Brown   PetscBool      scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE,btheap = PETSC_FALSE;
195c66b693SKris Buschelman 
205c66b693SKris Buschelman   PetscFunctionBegin;
2126be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
22152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
230298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,NULL);CHKERRQ(ierr);
240298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,NULL);CHKERRQ(ierr);
250298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,NULL);CHKERRQ(ierr);
260298fd71SBarry Smith     ierr = PetscOptionsBool("-matmatmult_btheap","Use btheap implementation of symbolic factorization C=A*B","",btheap,&btheap,NULL);CHKERRQ(ierr);
27d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
28*3ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
293c50cad2SHong Zhang     if (scalable_fast) {
303c50cad2SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
31aacf854cSBarry Smith     } else if (scalable) {
32aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
330ced3a2bSJed Brown     } else if (heap) {
340ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
358a07c6f1SJed Brown     } else if (btheap) {
368a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
3725023028SHong Zhang     } else {
3826be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3925023028SHong Zhang     }
40*3ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
4126be0446SHong Zhang   }
425c913ed7SHong Zhang 
43*3ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
4401e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
45*3ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
465c66b693SKris Buschelman   PetscFunctionReturn(0);
475c66b693SKris Buschelman }
481c24bd37SHong Zhang 
49e9e4536cSHong Zhang #undef __FUNCT__
50b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
51b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
52b561aa9dSHong Zhang {
53b561aa9dSHong Zhang   PetscErrorCode     ierr;
54b561aa9dSHong Zhang   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
55971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
56c1ba5575SJed Brown   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
57b561aa9dSHong Zhang   PetscReal          afill;
58fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
59fb908031SHong Zhang   PetscBT            lnkbt;
600298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
61b561aa9dSHong Zhang 
62b561aa9dSHong Zhang   PetscFunctionBegin;
63bd958071SHong Zhang   /* Get ci and cj */
64bd958071SHong Zhang   /*---------------*/
65fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
66fb908031SHong Zhang   /* free space for accumulating nonzero column info */
67fb908031SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
68fb908031SHong Zhang   ci[0] = 0;
69fb908031SHong Zhang 
70fb908031SHong Zhang   /* create and initialize a linked list */
71fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
72fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
73fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
74fb908031SHong Zhang 
75fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
76fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
772205254eSKarl Rupp 
78fb908031SHong Zhang   current_space = free_space;
79fb908031SHong Zhang 
80fb908031SHong Zhang   /* Determine ci and cj */
81fb908031SHong Zhang   for (i=0; i<am; i++) {
82fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
83fb908031SHong Zhang     aj   = a->j + ai[i];
84fb908031SHong Zhang     for (j=0; j<anzi; j++) {
85fb908031SHong Zhang       brow = aj[j];
86fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
87fb908031SHong Zhang       bj   = b->j + bi[brow];
88fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
89fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
90fb908031SHong Zhang     }
91fb908031SHong Zhang     cnzi = lnk[0];
92fb908031SHong Zhang 
93fb908031SHong Zhang     /* If free space is not available, make more free space */
94fb908031SHong Zhang     /* Double the amount of total space in the list */
95fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
96fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
97fb908031SHong Zhang       ndouble++;
98fb908031SHong Zhang     }
99fb908031SHong Zhang 
100fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
101fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1022205254eSKarl Rupp 
103fb908031SHong Zhang     current_space->array           += cnzi;
104fb908031SHong Zhang     current_space->local_used      += cnzi;
105fb908031SHong Zhang     current_space->local_remaining -= cnzi;
1062205254eSKarl Rupp 
107fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
108fb908031SHong Zhang   }
109fb908031SHong Zhang 
110fb908031SHong Zhang   /* Column indices are in the list of free space */
111fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
112fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
113fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
114fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
115fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
116b561aa9dSHong Zhang 
11726be0446SHong Zhang   /* put together the new symbolic matrix */
118ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1192205254eSKarl Rupp 
120a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
121a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
12258c24d83SHong Zhang 
12358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
12458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
12558c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
126aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
127e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
12858c24d83SHong Zhang   c->nonew                  = 0;
129002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1300b7e3e3dSHong Zhang 
1317212da7cSHong Zhang   /* set MatInfo */
1327212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
133f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1347212da7cSHong Zhang   c->maxnz                     = ci[am];
1357212da7cSHong Zhang   c->nz                        = ci[am];
136fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1377212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1387212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1397212da7cSHong Zhang 
1407212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1417212da7cSHong Zhang   if (ci[am]) {
142fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
143f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
144f2b054eeSHong Zhang   } else {
145f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
146be0fcf8dSHong Zhang   }
147f2b054eeSHong Zhang #endif
14858c24d83SHong Zhang   PetscFunctionReturn(0);
14958c24d83SHong Zhang }
150d50806bdSBarry Smith 
15126be0446SHong Zhang #undef __FUNCT__
15226be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
153dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
154d50806bdSBarry Smith {
155dfbe8321SBarry Smith   PetscErrorCode ierr;
156d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
157d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
158d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
159d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
16038baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
161c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
162519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
163aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
164aa1aec99SHong Zhang   PetscScalar    *ab_dense;
165d50806bdSBarry Smith 
166d50806bdSBarry Smith   PetscFunctionBegin;
167aa1aec99SHong Zhang   /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */
168aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
169aa1aec99SHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
170aa1aec99SHong Zhang     c->a      = ca;
171aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
172aa1aec99SHong Zhang 
173aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
174aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
1752205254eSKarl Rupp 
176aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
177aa1aec99SHong Zhang   } else {
178aa1aec99SHong Zhang     ca       = c->a;
179aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
180aa1aec99SHong Zhang   }
181aa1aec99SHong Zhang 
182c124e916SHong Zhang   /* clean old values in C */
183c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
184d50806bdSBarry Smith   /* Traverse A row-wise. */
185d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
186d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
187d50806bdSBarry Smith   for (i=0; i<am; i++) {
188d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
189d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
190519eb980SHong Zhang       brow = aj[j];
191d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
192d50806bdSBarry Smith       bjj  = bj + bi[brow];
193d50806bdSBarry Smith       baj  = ba + bi[brow];
194519eb980SHong Zhang       /* perform dense axpy */
19536ec6d2dSHong Zhang       valtmp = aa[j];
196519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
19736ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
198519eb980SHong Zhang       }
199519eb980SHong Zhang       flops += 2*bnzi;
200519eb980SHong Zhang     }
201c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
202c58ca2e3SHong Zhang 
203c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
204519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
205519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
206519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
207519eb980SHong Zhang     }
208519eb980SHong Zhang     flops += cnzi;
209519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
210519eb980SHong Zhang   }
211c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
212c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
213c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
214c58ca2e3SHong Zhang   PetscFunctionReturn(0);
215c58ca2e3SHong Zhang }
216c58ca2e3SHong Zhang 
217c58ca2e3SHong Zhang #undef __FUNCT__
21825023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
21925023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
220c58ca2e3SHong Zhang {
221c58ca2e3SHong Zhang   PetscErrorCode ierr;
222c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
223c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
224c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
225c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
226c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
227c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
228c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
22936ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
230c58ca2e3SHong Zhang   PetscInt       nextb;
231c58ca2e3SHong Zhang 
232c58ca2e3SHong Zhang   PetscFunctionBegin;
233c58ca2e3SHong Zhang   /* clean old values in C */
234c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
235c58ca2e3SHong Zhang   /* Traverse A row-wise. */
236c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
237c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
238519eb980SHong Zhang   for (i=0; i<am; i++) {
239519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
240519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
241519eb980SHong Zhang     for (j=0; j<anzi; j++) {
242519eb980SHong Zhang       brow = aj[j];
243519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
244519eb980SHong Zhang       bjj  = bj + bi[brow];
245519eb980SHong Zhang       baj  = ba + bi[brow];
246519eb980SHong Zhang       /* perform sparse axpy */
24736ec6d2dSHong Zhang       valtmp = aa[j];
248c124e916SHong Zhang       nextb  = 0;
249c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
250c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
25136ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
252c124e916SHong Zhang         }
253d50806bdSBarry Smith       }
254d50806bdSBarry Smith       flops += 2*bnzi;
255d50806bdSBarry Smith     }
256519eb980SHong Zhang     aj += anzi; aa += anzi;
257519eb980SHong Zhang     cj += cnzi; ca += cnzi;
258d50806bdSBarry Smith   }
259c58ca2e3SHong Zhang 
260716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
262d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
263d50806bdSBarry Smith   PetscFunctionReturn(0);
264d50806bdSBarry Smith }
265bc011b1eSHong Zhang 
266e9e4536cSHong Zhang #undef __FUNCT__
2673c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2683c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
26925296bd5SBarry Smith {
27025296bd5SBarry Smith   PetscErrorCode     ierr;
27125296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
27225296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
2733c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
27425296bd5SBarry Smith   MatScalar          *ca;
27525296bd5SBarry Smith   PetscReal          afill;
2763c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2770298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
27825296bd5SBarry Smith 
27925296bd5SBarry Smith   PetscFunctionBegin;
2803c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2813c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2823c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2833c50cad2SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2843c50cad2SHong Zhang   ci[0] = 0;
2853c50cad2SHong Zhang 
2863c50cad2SHong Zhang   /* create and initialize a linked list */
2873c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2883c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2893c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
2903c50cad2SHong Zhang 
2913c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
2923c50cad2SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
2933c50cad2SHong Zhang   current_space = free_space;
2943c50cad2SHong Zhang 
2953c50cad2SHong Zhang   /* Determine ci and cj */
2963c50cad2SHong Zhang   for (i=0; i<am; i++) {
2973c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
2983c50cad2SHong Zhang     aj   = a->j + ai[i];
2993c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3003c50cad2SHong Zhang       brow = aj[j];
3013c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3023c50cad2SHong Zhang       bj   = b->j + bi[brow];
3033c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3043c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3053c50cad2SHong Zhang     }
3063c50cad2SHong Zhang     cnzi = lnk[1];
3073c50cad2SHong Zhang 
3083c50cad2SHong Zhang     /* If free space is not available, make more free space */
3093c50cad2SHong Zhang     /* Double the amount of total space in the list */
3103c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3113c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3123c50cad2SHong Zhang       ndouble++;
3133c50cad2SHong Zhang     }
3143c50cad2SHong Zhang 
3153c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3163c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3172205254eSKarl Rupp 
3183c50cad2SHong Zhang     current_space->array           += cnzi;
3193c50cad2SHong Zhang     current_space->local_used      += cnzi;
3203c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3212205254eSKarl Rupp 
3223c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3233c50cad2SHong Zhang   }
3243c50cad2SHong Zhang 
3253c50cad2SHong Zhang   /* Column indices are in the list of free space */
3263c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3273c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3283c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3293c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3303c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
33125296bd5SBarry Smith 
33225296bd5SBarry Smith   /* Allocate space for ca */
33325296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
33425296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
33525296bd5SBarry Smith 
33625296bd5SBarry Smith   /* put together the new symbolic matrix */
337ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
3382205254eSKarl Rupp 
339a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
340a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
34125296bd5SBarry Smith 
34225296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
34325296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
34425296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
34525296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
34625296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
34725296bd5SBarry Smith   c->nonew   = 0;
3482205254eSKarl Rupp 
34925296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
35025296bd5SBarry Smith 
35125296bd5SBarry Smith   /* set MatInfo */
35225296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
35325296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
35425296bd5SBarry Smith   c->maxnz                     = ci[am];
35525296bd5SBarry Smith   c->nz                        = ci[am];
3563c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
35725296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
35825296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
35925296bd5SBarry Smith 
36025296bd5SBarry Smith #if defined(PETSC_USE_INFO)
36125296bd5SBarry Smith   if (ci[am]) {
3623c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
36325296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
36425296bd5SBarry Smith   } else {
36525296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
36625296bd5SBarry Smith   }
36725296bd5SBarry Smith #endif
36825296bd5SBarry Smith   PetscFunctionReturn(0);
36925296bd5SBarry Smith }
37025296bd5SBarry Smith 
37125296bd5SBarry Smith 
37225296bd5SBarry Smith #undef __FUNCT__
37325023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
37425023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
375e9e4536cSHong Zhang {
376e9e4536cSHong Zhang   PetscErrorCode     ierr;
377e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
378bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
37925c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
380e9e4536cSHong Zhang   MatScalar          *ca;
381e9e4536cSHong Zhang   PetscReal          afill;
382bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
3830298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
384e9e4536cSHong Zhang 
385e9e4536cSHong Zhang   PetscFunctionBegin;
386bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
387bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
388bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
389bd958071SHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
390bd958071SHong Zhang   ci[0] = 0;
391bd958071SHong Zhang 
392bd958071SHong Zhang   /* create and initialize a linked list */
393bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
394bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
395bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
396bd958071SHong Zhang 
397bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
398bd958071SHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
399bd958071SHong Zhang   current_space = free_space;
400bd958071SHong Zhang 
401bd958071SHong Zhang   /* Determine ci and cj */
402bd958071SHong Zhang   for (i=0; i<am; i++) {
403bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
404bd958071SHong Zhang     aj   = a->j + ai[i];
405bd958071SHong Zhang     for (j=0; j<anzi; j++) {
406bd958071SHong Zhang       brow = aj[j];
407bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
408bd958071SHong Zhang       bj   = b->j + bi[brow];
409bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
410bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
411bd958071SHong Zhang     }
412bd958071SHong Zhang     cnzi = lnk[0];
413bd958071SHong Zhang 
414bd958071SHong Zhang     /* If free space is not available, make more free space */
415bd958071SHong Zhang     /* Double the amount of total space in the list */
416bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
417bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
418bd958071SHong Zhang       ndouble++;
419bd958071SHong Zhang     }
420bd958071SHong Zhang 
421bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
422bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4232205254eSKarl Rupp 
424bd958071SHong Zhang     current_space->array           += cnzi;
425bd958071SHong Zhang     current_space->local_used      += cnzi;
426bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4272205254eSKarl Rupp 
428bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
429bd958071SHong Zhang   }
430bd958071SHong Zhang 
431bd958071SHong Zhang   /* Column indices are in the list of free space */
432bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
433bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
434bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
435bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
436bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
437e9e4536cSHong Zhang 
438e9e4536cSHong Zhang   /* Allocate space for ca */
439bd958071SHong Zhang   /*-----------------------*/
440e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
441e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
442e9e4536cSHong Zhang 
443e9e4536cSHong Zhang   /* put together the new symbolic matrix */
444ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
4452205254eSKarl Rupp 
446a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
447a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
448e9e4536cSHong Zhang 
449e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
450e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
451e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
452e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
453e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
454e9e4536cSHong Zhang   c->nonew   = 0;
4552205254eSKarl Rupp 
45625023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
457e9e4536cSHong Zhang 
458e9e4536cSHong Zhang   /* set MatInfo */
459e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
460e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
461e9e4536cSHong Zhang   c->maxnz                     = ci[am];
462e9e4536cSHong Zhang   c->nz                        = ci[am];
463bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
464e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
465e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
466e9e4536cSHong Zhang 
467e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
468e9e4536cSHong Zhang   if (ci[am]) {
469bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
470e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
471e9e4536cSHong Zhang   } else {
472e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
473e9e4536cSHong Zhang   }
474e9e4536cSHong Zhang #endif
475e9e4536cSHong Zhang   PetscFunctionReturn(0);
476e9e4536cSHong Zhang }
477e9e4536cSHong Zhang 
4780ced3a2bSJed Brown #undef __FUNCT__
4790ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
4800ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
4810ced3a2bSJed Brown {
4820ced3a2bSJed Brown   PetscErrorCode     ierr;
4830ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4840ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
4850ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
4860ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
4870ced3a2bSJed Brown   PetscReal          afill;
4880ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
4890298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
4900ced3a2bSJed Brown   PetscHeap          h;
4910ced3a2bSJed Brown 
4920ced3a2bSJed Brown   PetscFunctionBegin;
4930ced3a2bSJed Brown   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
4940ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
4950ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
4960ced3a2bSJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4970ced3a2bSJed Brown   ci[0] = 0;
4980ced3a2bSJed Brown 
4990ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
5000ced3a2bSJed Brown   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5010ced3a2bSJed Brown   current_space = free_space;
5020ced3a2bSJed Brown 
5030ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
5040ced3a2bSJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
5050ced3a2bSJed Brown 
5060ced3a2bSJed Brown   /* Determine ci and cj */
5070ced3a2bSJed Brown   for (i=0; i<am; i++) {
5080ced3a2bSJed 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 */
5090ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5100ced3a2bSJed Brown     ci[i+1] = ci[i];
5110ced3a2bSJed Brown     /* Populate the min heap */
5120ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5130ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5140ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5150ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5160ced3a2bSJed Brown       }
5170ced3a2bSJed Brown     }
5180ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5190ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5200ced3a2bSJed Brown     while (j >= 0) {
5210ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
5220ced3a2bSJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
5230ced3a2bSJed Brown         ndouble++;
5240ced3a2bSJed Brown       }
5250ced3a2bSJed Brown       *(current_space->array++) = col;
5260ced3a2bSJed Brown       current_space->local_used++;
5270ced3a2bSJed Brown       current_space->local_remaining--;
5280ced3a2bSJed Brown       ci[i+1]++;
5290ced3a2bSJed Brown 
5300ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5310ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5320ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5330ced3a2bSJed Brown         PetscInt j2,col2;
5340ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5350ced3a2bSJed Brown         if (col2 != col) break;
5360ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5370ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5380ced3a2bSJed Brown       }
5390ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5400ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5410ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5420ced3a2bSJed Brown     }
5430ced3a2bSJed Brown   }
5440ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5450ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5460ced3a2bSJed Brown 
5470ced3a2bSJed Brown   /* Column indices are in the list of free space */
5480ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5490ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
5500ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5510ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5520ced3a2bSJed Brown 
5530ced3a2bSJed Brown   /* put together the new symbolic matrix */
554ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
5552205254eSKarl Rupp 
5560ced3a2bSJed Brown   (*C)->rmap->bs = A->rmap->bs;
5570ced3a2bSJed Brown   (*C)->cmap->bs = B->cmap->bs;
5580ced3a2bSJed Brown 
5590ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5600ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5610ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5620ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5630ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5640ced3a2bSJed Brown   c->nonew   = 0;
56526fbe8dcSKarl Rupp 
56689d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5670ced3a2bSJed Brown 
5680ced3a2bSJed Brown   /* set MatInfo */
5690ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5700ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5710ced3a2bSJed Brown   c->maxnz                     = ci[am];
5720ced3a2bSJed Brown   c->nz                        = ci[am];
5730ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5740ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5750ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5760ced3a2bSJed Brown 
5770ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5780ced3a2bSJed Brown   if (ci[am]) {
5790ced3a2bSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
5800ced3a2bSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
5810ced3a2bSJed Brown   } else {
5820ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
5830ced3a2bSJed Brown   }
5840ced3a2bSJed Brown #endif
5850ced3a2bSJed Brown   PetscFunctionReturn(0);
5860ced3a2bSJed Brown }
587e9e4536cSHong Zhang 
5888a07c6f1SJed Brown #undef __FUNCT__
5898a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
5908a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
5918a07c6f1SJed Brown {
5928a07c6f1SJed Brown   PetscErrorCode     ierr;
5938a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5948a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
5958a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
5968a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
5978a07c6f1SJed Brown   PetscReal          afill;
5988a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
5990298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6008a07c6f1SJed Brown   PetscHeap          h;
6018a07c6f1SJed Brown   PetscBT            bt;
6028a07c6f1SJed Brown 
6038a07c6f1SJed Brown   PetscFunctionBegin;
6048a07c6f1SJed Brown   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
6058a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6068a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
6078a07c6f1SJed Brown   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6088a07c6f1SJed Brown   ci[0] = 0;
6098a07c6f1SJed Brown 
6108a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
6118a07c6f1SJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6122205254eSKarl Rupp 
6138a07c6f1SJed Brown   current_space = free_space;
6148a07c6f1SJed Brown 
6158a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
6168a07c6f1SJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
6178a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6188a07c6f1SJed Brown 
6198a07c6f1SJed Brown   /* Determine ci and cj */
6208a07c6f1SJed Brown   for (i=0; i<am; i++) {
6218a07c6f1SJed 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 */
6228a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6238a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6248a07c6f1SJed Brown     ci[i+1] = ci[i];
6258a07c6f1SJed Brown     /* Populate the min heap */
6268a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6278a07c6f1SJed Brown       PetscInt brow = acol[j];
6288a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6298a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6308a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6318a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6328a07c6f1SJed Brown           bb[j]++;
6338a07c6f1SJed Brown           break;
6348a07c6f1SJed Brown         }
6358a07c6f1SJed Brown       }
6368a07c6f1SJed Brown     }
6378a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6388a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6398a07c6f1SJed Brown     while (j >= 0) {
6408a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6410298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
6428a07c6f1SJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
6438a07c6f1SJed Brown         ndouble++;
6448a07c6f1SJed Brown       }
6458a07c6f1SJed Brown       *(current_space->array++) = col;
6468a07c6f1SJed Brown       current_space->local_used++;
6478a07c6f1SJed Brown       current_space->local_remaining--;
6488a07c6f1SJed Brown       ci[i+1]++;
6498a07c6f1SJed Brown 
6508a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6518a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6528a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6538a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6548a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6558a07c6f1SJed Brown           bb[j]++;
6568a07c6f1SJed Brown           break;
6578a07c6f1SJed Brown         }
6588a07c6f1SJed Brown       }
6598a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6608a07c6f1SJed Brown     }
6618a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6628a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6638a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6648a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6658a07c6f1SJed Brown     }
6668a07c6f1SJed Brown   }
6678a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6688a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6698a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6708a07c6f1SJed Brown 
6718a07c6f1SJed Brown   /* Column indices are in the list of free space */
6728a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6738a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
6748a07c6f1SJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6758a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6768a07c6f1SJed Brown 
6778a07c6f1SJed Brown   /* put together the new symbolic matrix */
678ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
6792205254eSKarl Rupp 
6808a07c6f1SJed Brown   (*C)->rmap->bs = A->rmap->bs;
6818a07c6f1SJed Brown   (*C)->cmap->bs = B->cmap->bs;
6828a07c6f1SJed Brown 
6838a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6848a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6858a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
6868a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
6878a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
6888a07c6f1SJed Brown   c->nonew   = 0;
68926fbe8dcSKarl Rupp 
69089d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
6918a07c6f1SJed Brown 
6928a07c6f1SJed Brown   /* set MatInfo */
6938a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6948a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
6958a07c6f1SJed Brown   c->maxnz                     = ci[am];
6968a07c6f1SJed Brown   c->nz                        = ci[am];
6978a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
6988a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
6998a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7008a07c6f1SJed Brown 
7018a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7028a07c6f1SJed Brown   if (ci[am]) {
7038a07c6f1SJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
7048a07c6f1SJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
7058a07c6f1SJed Brown   } else {
7068a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7078a07c6f1SJed Brown   }
7088a07c6f1SJed Brown #endif
7098a07c6f1SJed Brown   PetscFunctionReturn(0);
7108a07c6f1SJed Brown }
7118a07c6f1SJed Brown 
712d2853540SHong Zhang /* This routine is not used. Should be removed! */
713bc011b1eSHong Zhang #undef __FUNCT__
7146fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
7156fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
7165df89d91SHong Zhang {
717bc011b1eSHong Zhang   PetscErrorCode ierr;
718bc011b1eSHong Zhang 
719bc011b1eSHong Zhang   PetscFunctionBegin;
720bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
721*3ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
7226fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
723*3ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
724bc011b1eSHong Zhang   }
725*3ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
7266fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
727*3ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
728bc011b1eSHong Zhang   PetscFunctionReturn(0);
729bc011b1eSHong Zhang }
730bc011b1eSHong Zhang 
731bc011b1eSHong Zhang #undef __FUNCT__
732e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
733e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
7342128a86cSHong Zhang {
7352128a86cSHong Zhang   PetscErrorCode      ierr;
736e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
7372128a86cSHong Zhang 
7382128a86cSHong Zhang   PetscFunctionBegin;
739b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
7402128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
7412128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
7422128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
7432128a86cSHong Zhang   PetscFunctionReturn(0);
7442128a86cSHong Zhang }
7452128a86cSHong Zhang 
7462128a86cSHong Zhang #undef __FUNCT__
7472128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
7482128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
7492128a86cSHong Zhang {
7502128a86cSHong Zhang   PetscErrorCode      ierr;
7512128a86cSHong Zhang   PetscContainer      container;
7520298fd71SBarry Smith   Mat_MatMatTransMult *multtrans=NULL;
7532128a86cSHong Zhang 
7542128a86cSHong Zhang   PetscFunctionBegin;
755e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr);
7562128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
7572128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr);
7582205254eSKarl Rupp 
7592128a86cSHong Zhang   A->ops->destroy = multtrans->destroy;
7602128a86cSHong Zhang   if (A->ops->destroy) {
7612128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
7622128a86cSHong Zhang   }
763e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
7642128a86cSHong Zhang   PetscFunctionReturn(0);
7652128a86cSHong Zhang }
7662128a86cSHong Zhang 
7672128a86cSHong Zhang #undef __FUNCT__
7686fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
7696fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
770bc011b1eSHong Zhang {
7715df89d91SHong Zhang   PetscErrorCode      ierr;
77281d82fe4SHong Zhang   Mat                 Bt;
77381d82fe4SHong Zhang   PetscInt            *bti,*btj;
774e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
7752128a86cSHong Zhang   PetscContainer      container;
776d2853540SHong Zhang 
77781d82fe4SHong Zhang   PetscFunctionBegin;
77881d82fe4SHong Zhang   /* create symbolic Bt */
77981d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
7800298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
7812205254eSKarl Rupp 
782c3e5699bSSatish Balay   Bt->rmap->bs = A->cmap->bs;
783c3e5699bSSatish Balay   Bt->cmap->bs = B->cmap->bs;
78481d82fe4SHong Zhang 
78581d82fe4SHong Zhang   /* get symbolic C=A*Bt */
78681d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
78781d82fe4SHong Zhang 
7882128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
789e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
7902128a86cSHong Zhang 
7912128a86cSHong Zhang   /* attach the supporting struct to C */
7922128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
7932128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
794e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
795e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
7962128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
7972128a86cSHong Zhang 
7982128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
7992128a86cSHong Zhang   multtrans->destroy     = (*C)->ops->destroy;
8002128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8012128a86cSHong Zhang 
8020298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-matmattransmult_color",&multtrans->usecoloring,NULL);CHKERRQ(ierr);
8032128a86cSHong Zhang   if (multtrans->usecoloring) {
804b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
805b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
8062128a86cSHong Zhang     ISColoring           iscoloring;
8072128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
808e8354b3eSHong Zhang 
809*3ff4c91cSHong Zhang     //ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
810*3ff4c91cSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr);
811b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
8122205254eSKarl Rupp 
8132128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
8142205254eSKarl Rupp 
8152128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
8162128a86cSHong Zhang 
8172128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
8182128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
8192128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
8202128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
8210298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
8222205254eSKarl Rupp 
8232128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
8242128a86cSHong Zhang     multtrans->Bt_den   = Bt_dense;
8252128a86cSHong Zhang 
8262128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
8272128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
8282128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
8290298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
8302205254eSKarl Rupp 
8312128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
8322128a86cSHong Zhang     multtrans->ABt_den  = C_dense;
833f94ccd6cSHong Zhang 
834f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
8351ce71dffSSatish Balay     {
836f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
83722d28d08SBarry Smith       ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr);
8381ce71dffSSatish Balay     }
839f94ccd6cSHong Zhang #endif
8402128a86cSHong Zhang   }
841e99be685SHong Zhang   /* clean up */
842e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
843e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8442128a86cSHong Zhang 
845f94ccd6cSHong Zhang 
846f94ccd6cSHong Zhang 
84781d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
84881d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
8490298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
8501fa1209cSHong Zhang   Mat_SeqAIJ         *a        =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
8511fa1209cSHong Zhang   PetscInt           *ai       =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
8521fa1209cSHong Zhang   PetscInt           am        =A->rmap->N,bm=B->rmap->N;
8531fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
8541fa1209cSHong Zhang   MatScalar          *ca;
8551fa1209cSHong Zhang   PetscBT            lnkbt;
8561fa1209cSHong Zhang   PetscReal          afill;
8575df89d91SHong Zhang 
8581fa1209cSHong Zhang   /* Allocate row pointer array ci  */
8591fa1209cSHong Zhang   ierr  = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
8601fa1209cSHong Zhang   ci[0] = 0;
8611fa1209cSHong Zhang 
8621fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
8631fa1209cSHong Zhang   nlnk = bm+1;
8641fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
8651fa1209cSHong Zhang 
8661fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
8671fa1209cSHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
8681fa1209cSHong Zhang   current_space = free_space;
8691fa1209cSHong Zhang 
8701fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
8711fa1209cSHong Zhang   for (i=0; i<am; i++) {
8721fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
8731fa1209cSHong Zhang     cnzi = 0;
8741fa1209cSHong Zhang     acol = aj + ai[i];
8751fa1209cSHong Zhang     for (j=0; j<bm; j++) {
8761fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
8771fa1209cSHong Zhang       bcol = bj + bi[j];
87881d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
8791fa1209cSHong Zhang       ka = 0; kb = 0;
8801fa1209cSHong Zhang       while (ka < anzi && kb < bnzj) {
88181d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
88281d82fe4SHong Zhang         if (ka == anzi) break;
88381d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
88481d82fe4SHong Zhang         if (kb == bnzj) break;
8851fa1209cSHong Zhang         if (acol[ka] == bcol[kb]) { /* add nonzero c(i,j) to lnk */
8861fa1209cSHong Zhang           index[0] = j;
8871fa1209cSHong Zhang           ierr     = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
8881fa1209cSHong Zhang           cnzi++;
8891fa1209cSHong Zhang           break;
8901fa1209cSHong Zhang         }
8911fa1209cSHong Zhang       }
8921fa1209cSHong Zhang     }
8931fa1209cSHong Zhang 
8941fa1209cSHong Zhang     /* If free space is not available, make more free space */
8951fa1209cSHong Zhang     /* Double the amount of total space in the list */
8961fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
8971fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
8981fa1209cSHong Zhang       nspacedouble++;
8991fa1209cSHong Zhang     }
9001fa1209cSHong Zhang 
9011fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
9021fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
9032205254eSKarl Rupp 
9041fa1209cSHong Zhang     current_space->array           += cnzi;
9051fa1209cSHong Zhang     current_space->local_used      += cnzi;
9061fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
9071fa1209cSHong Zhang 
9081fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
9091fa1209cSHong Zhang   }
9101fa1209cSHong Zhang 
9111fa1209cSHong Zhang 
9121fa1209cSHong Zhang   /* Column indices are in the list of free space.
9131fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
9141fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
9151fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
9161fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
9171fa1209cSHong Zhang 
9181fa1209cSHong Zhang   /* Allocate space for ca */
9191fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
9201fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
9211fa1209cSHong Zhang 
9221fa1209cSHong Zhang   /* put together the new symbolic matrix */
923ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bm,ci,cj,ca,C);CHKERRQ(ierr);
9242205254eSKarl Rupp 
925a2f3521dSMark F. Adams   (*C)->rmap->bs = A->cmap->bs;
926a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
9271fa1209cSHong Zhang 
9281fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
9291fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
9301fa1209cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
9311fa1209cSHong Zhang   c->free_a  = PETSC_TRUE;
9321fa1209cSHong Zhang   c->free_ij = PETSC_TRUE;
9331fa1209cSHong Zhang   c->nonew   = 0;
9341fa1209cSHong Zhang 
9351fa1209cSHong Zhang   /* set MatInfo */
9361fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
9371fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
9381fa1209cSHong Zhang   c->maxnz                     = ci[am];
9391fa1209cSHong Zhang   c->nz                        = ci[am];
9401fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
9411fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
9421fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
9431fa1209cSHong Zhang 
9441fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
9451fa1209cSHong Zhang   if (ci[am]) {
9461fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
9476fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
9481fa1209cSHong Zhang   } else {
9491fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
9501fa1209cSHong Zhang   }
9511fa1209cSHong Zhang #endif
95281d82fe4SHong Zhang #endif
9535df89d91SHong Zhang   PetscFunctionReturn(0);
9545df89d91SHong Zhang }
9555df89d91SHong Zhang 
9566973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
9575df89d91SHong Zhang #undef __FUNCT__
9586fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9596fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9605df89d91SHong Zhang {
9615df89d91SHong Zhang   PetscErrorCode      ierr;
9625df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
963e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
96481d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9655df89d91SHong Zhang   PetscLogDouble      flops=0.0;
966aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
967e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
9682128a86cSHong Zhang   PetscContainer      container;
9696973769bSHong Zhang #if defined(USE_ARRAY)
9706973769bSHong Zhang   MatScalar *spdot;
9716973769bSHong Zhang #endif
9725df89d91SHong Zhang 
9735df89d91SHong Zhang   PetscFunctionBegin;
97458ed3ceaSHong Zhang   /* clear old values in C */
97558ed3ceaSHong Zhang   if (!c->a) {
97658ed3ceaSHong Zhang     ierr      = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
97758ed3ceaSHong Zhang     c->a      = ca;
97858ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
97958ed3ceaSHong Zhang   } else {
98058ed3ceaSHong Zhang     ca =  c->a;
98158ed3ceaSHong Zhang   }
98258ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
98358ed3ceaSHong Zhang 
984e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr);
9852128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
9862128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr);
9872128a86cSHong Zhang   if (multtrans->usecoloring) {
988b9af6bddSHong Zhang     MatTransposeColoring matcoloring = multtrans->matcoloring;
989c8db22f6SHong Zhang     Mat                  Bt_dense;
990c8db22f6SHong Zhang     PetscInt             m,n;
991a0b95ffbSSatish Balay     Mat                  C_dense = multtrans->ABt_den;
992a0b95ffbSSatish Balay 
9932128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
994c8db22f6SHong Zhang     ierr     = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
995c8db22f6SHong Zhang 
996b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
997b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
998c8db22f6SHong Zhang 
999c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
10002128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
1001c8db22f6SHong Zhang 
10022128a86cSHong Zhang     /* Recover C from C_dense */
1003b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
1004c8db22f6SHong Zhang     PetscFunctionReturn(0);
1005c8db22f6SHong Zhang   }
1006c8db22f6SHong Zhang 
10076973769bSHong Zhang #if defined(USE_ARRAY)
10086973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
1009e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
1010e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
10116973769bSHong Zhang #endif
10126973769bSHong Zhang 
10131fa1209cSHong Zhang   for (i=0; i<cm; i++) {
101481d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
101581d82fe4SHong Zhang     acol = aj + ai[i];
10166973769bSHong Zhang     aval = aa + ai[i];
10171fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
10181fa1209cSHong Zhang     ccol = cj + ci[i];
10196973769bSHong Zhang     cval = ca + ci[i];
10201fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
102181d82fe4SHong Zhang       brow = ccol[j];
102281d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
102381d82fe4SHong Zhang       bcol = bj + bi[brow];
10246973769bSHong Zhang       bval = ba + bi[brow];
10256973769bSHong Zhang 
102681d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
10276973769bSHong Zhang #if defined(USE_ARRAY)
10286973769bSHong Zhang       /* put ba to spdot array */
10296973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
10306973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
10316973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++) {
10326973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
10336973769bSHong Zhang       }
10346973769bSHong Zhang       /* zero spdot array */
10356973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
10366973769bSHong Zhang #else
103781d82fe4SHong Zhang       nexta = 0; nextb = 0;
103881d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
10397b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
104081d82fe4SHong Zhang         if (nexta == anzi) break;
10417b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
104281d82fe4SHong Zhang         if (nextb == bnzj) break;
104381d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
10446973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
104581d82fe4SHong Zhang           nexta++; nextb++;
104681d82fe4SHong Zhang           flops += 2;
10471fa1209cSHong Zhang         }
10481fa1209cSHong Zhang       }
10496973769bSHong Zhang #endif
105081d82fe4SHong Zhang     }
105181d82fe4SHong Zhang   }
105281d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105381d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105481d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10556973769bSHong Zhang #if defined(USE_ARRAY)
10566973769bSHong Zhang   ierr = PetscFree(spdot);
10576973769bSHong Zhang #endif
10585df89d91SHong Zhang   PetscFunctionReturn(0);
10595df89d91SHong Zhang }
10605df89d91SHong Zhang 
10615df89d91SHong Zhang #undef __FUNCT__
106275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
10630adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10640adebc6cSBarry Smith {
10655df89d91SHong Zhang   PetscErrorCode ierr;
10665df89d91SHong Zhang 
10675df89d91SHong Zhang   PetscFunctionBegin;
10685df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
106975648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
10705df89d91SHong Zhang   }
107175648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
10725df89d91SHong Zhang   PetscFunctionReturn(0);
10735df89d91SHong Zhang }
10745df89d91SHong Zhang 
10755df89d91SHong Zhang #undef __FUNCT__
107675648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
107775648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10785df89d91SHong Zhang {
1079bc011b1eSHong Zhang   PetscErrorCode ierr;
1080bc011b1eSHong Zhang   Mat            At;
108138baddfdSBarry Smith   PetscInt       *ati,*atj;
1082bc011b1eSHong Zhang 
1083bc011b1eSHong Zhang   PetscFunctionBegin;
1084bc011b1eSHong Zhang   /* create symbolic At */
1085bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10860298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
10872205254eSKarl Rupp 
1088a2f3521dSMark F. Adams   At->rmap->bs = A->cmap->bs;
1089a2f3521dSMark F. Adams   At->cmap->bs = B->cmap->bs;
1090bc011b1eSHong Zhang 
1091bc011b1eSHong Zhang   /* get symbolic C=At*B */
1092bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1093bc011b1eSHong Zhang 
1094bc011b1eSHong Zhang   /* clean up */
10956bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1096bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1097bc011b1eSHong Zhang   PetscFunctionReturn(0);
1098bc011b1eSHong Zhang }
1099bc011b1eSHong Zhang 
1100bc011b1eSHong Zhang #undef __FUNCT__
110175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
110275648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1103bc011b1eSHong Zhang {
1104bc011b1eSHong Zhang   PetscErrorCode ierr;
11050fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1106d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1107d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1108d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1109aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1110bc011b1eSHong Zhang 
1111bc011b1eSHong Zhang   PetscFunctionBegin;
1112aa1aec99SHong Zhang   if (!c->a) {
1113aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
11142205254eSKarl Rupp 
1115aa1aec99SHong Zhang     c->a      = ca;
1116aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1117aa1aec99SHong Zhang   } else {
1118aa1aec99SHong Zhang     ca = c->a;
1119aa1aec99SHong Zhang   }
1120bc011b1eSHong Zhang   /* clear old values in C */
1121bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1122bc011b1eSHong Zhang 
1123bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1124bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1125bc011b1eSHong Zhang     bj   = b->j + bi[i];
1126bc011b1eSHong Zhang     ba   = b->a + bi[i];
1127bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1128bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1129bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1130bc011b1eSHong Zhang       nextb = 0;
11310fbc74f4SHong Zhang       crow  = *aj++;
1132bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1133bc011b1eSHong Zhang       caj   = ca + ci[crow];
1134bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1135bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
11360fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
11370fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1138bc011b1eSHong Zhang           nextb++;
1139bc011b1eSHong Zhang         }
1140bc011b1eSHong Zhang       }
1141bc011b1eSHong Zhang       flops += 2*bnzi;
11420fbc74f4SHong Zhang       aa++;
1143bc011b1eSHong Zhang     }
1144bc011b1eSHong Zhang   }
1145bc011b1eSHong Zhang 
1146bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1147bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1148bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1149bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1150bc011b1eSHong Zhang   PetscFunctionReturn(0);
1151bc011b1eSHong Zhang }
11529123193aSHong Zhang 
11539123193aSHong Zhang #undef __FUNCT__
11549123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
11559123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11569123193aSHong Zhang {
11579123193aSHong Zhang   PetscErrorCode ierr;
11589123193aSHong Zhang 
11599123193aSHong Zhang   PetscFunctionBegin;
11609123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
1161*3ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11629123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1163*3ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11649123193aSHong Zhang   }
1165*3ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11669123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
1167*3ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
11689123193aSHong Zhang   PetscFunctionReturn(0);
11699123193aSHong Zhang }
1170edd81eecSMatthew Knepley 
11719123193aSHong Zhang #undef __FUNCT__
11729123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11739123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11749123193aSHong Zhang {
11759123193aSHong Zhang   PetscErrorCode ierr;
11769123193aSHong Zhang 
11779123193aSHong Zhang   PetscFunctionBegin;
11785a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11792205254eSKarl Rupp 
1180d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11819123193aSHong Zhang   PetscFunctionReturn(0);
11829123193aSHong Zhang }
11839123193aSHong Zhang 
11849123193aSHong Zhang #undef __FUNCT__
11859123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11869123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11879123193aSHong Zhang {
1188f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11899123193aSHong Zhang   PetscErrorCode ierr;
1190dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1191dd6ea824SBarry Smith   MatScalar      *aa;
1192d0f46423SBarry Smith   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
1193f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
11949123193aSHong Zhang 
11959123193aSHong Zhang   PetscFunctionBegin;
1196f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1197e32f2f54SBarry 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);
1198e32f2f54SBarry 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);
1199e32f2f54SBarry 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);
12008c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12018c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1202f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1203f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1204f32d5d43SBarry Smith     colam = col*am;
1205f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1206f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1207f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1208f32d5d43SBarry Smith       aj = a->j + a->i[i];
1209f32d5d43SBarry Smith       aa = a->a + a->i[i];
1210f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1211f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
1212f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
1213f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
1214f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
12159123193aSHong Zhang       }
1216f32d5d43SBarry Smith       c[colam + i]       = r1;
1217f32d5d43SBarry Smith       c[colam + am + i]  = r2;
1218f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
1219f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
1220f32d5d43SBarry Smith     }
1221f32d5d43SBarry Smith     b1 += bm4;
1222f32d5d43SBarry Smith     b2 += bm4;
1223f32d5d43SBarry Smith     b3 += bm4;
1224f32d5d43SBarry Smith     b4 += bm4;
1225f32d5d43SBarry Smith   }
1226f32d5d43SBarry Smith   for (; col<cn; col++) {     /* over extra columns of C */
1227f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1228f32d5d43SBarry Smith       r1 = 0.0;
1229f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1230f32d5d43SBarry Smith       aj = a->j + a->i[i];
1231f32d5d43SBarry Smith       aa = a->a + a->i[i];
1232f32d5d43SBarry Smith 
1233f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1234f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
1235f32d5d43SBarry Smith       }
1236f32d5d43SBarry Smith       c[col*am + i] = r1;
1237f32d5d43SBarry Smith     }
1238f32d5d43SBarry Smith     b1 += bm;
1239f32d5d43SBarry Smith   }
1240dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
12418c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
12428c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1243f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1244f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1245f32d5d43SBarry Smith   PetscFunctionReturn(0);
1246f32d5d43SBarry Smith }
1247f32d5d43SBarry Smith 
1248f32d5d43SBarry Smith /*
12494324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1250f32d5d43SBarry Smith */
1251f32d5d43SBarry Smith #undef __FUNCT__
1252f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1253f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1254f32d5d43SBarry Smith {
1255f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1256f32d5d43SBarry Smith   PetscErrorCode ierr;
1257dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1258dd6ea824SBarry Smith   MatScalar      *aa;
1259d0f46423SBarry 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;
12604324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1261f32d5d43SBarry Smith 
1262f32d5d43SBarry Smith   PetscFunctionBegin;
1263f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
12648c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12658c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1266f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12674324174eSBarry Smith 
12684324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12694324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12704324174eSBarry Smith       colam = col*am;
12714324174eSBarry Smith       arm   = a->compressedrow.nrows;
12724324174eSBarry Smith       ii    = a->compressedrow.i;
12734324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12744324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12754324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12764324174eSBarry Smith         n  = ii[i+1] - ii[i];
12774324174eSBarry Smith         aj = a->j + ii[i];
12784324174eSBarry Smith         aa = a->a + ii[i];
12794324174eSBarry Smith         for (j=0; j<n; j++) {
12804324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12814324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12824324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12834324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12844324174eSBarry Smith         }
12854324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12864324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12874324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12884324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12894324174eSBarry Smith       }
12904324174eSBarry Smith       b1 += bm4;
12914324174eSBarry Smith       b2 += bm4;
12924324174eSBarry Smith       b3 += bm4;
12934324174eSBarry Smith       b4 += bm4;
12944324174eSBarry Smith     }
12954324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12964324174eSBarry Smith       colam = col*am;
12974324174eSBarry Smith       arm   = a->compressedrow.nrows;
12984324174eSBarry Smith       ii    = a->compressedrow.i;
12994324174eSBarry Smith       ridx  = a->compressedrow.rindex;
13004324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
13014324174eSBarry Smith         r1 = 0.0;
13024324174eSBarry Smith         n  = ii[i+1] - ii[i];
13034324174eSBarry Smith         aj = a->j + ii[i];
13044324174eSBarry Smith         aa = a->a + ii[i];
13054324174eSBarry Smith 
13064324174eSBarry Smith         for (j=0; j<n; j++) {
13074324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
13084324174eSBarry Smith         }
1309a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
13104324174eSBarry Smith       }
13114324174eSBarry Smith       b1 += bm;
13124324174eSBarry Smith     }
13134324174eSBarry Smith   } else {
1314f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1315f32d5d43SBarry Smith       colam = col*am;
1316f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1317f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1318f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1319f32d5d43SBarry Smith         aj = a->j + a->i[i];
1320f32d5d43SBarry Smith         aa = a->a + a->i[i];
1321f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1322f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1323f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1324f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1325f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1326f32d5d43SBarry Smith         }
1327f32d5d43SBarry Smith         c[colam + i]       += r1;
1328f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1329f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1330f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1331f32d5d43SBarry Smith       }
1332f32d5d43SBarry Smith       b1 += bm4;
1333f32d5d43SBarry Smith       b2 += bm4;
1334f32d5d43SBarry Smith       b3 += bm4;
1335f32d5d43SBarry Smith       b4 += bm4;
1336f32d5d43SBarry Smith     }
1337f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1338a2ea699eSBarry Smith       colam = col*am;
1339f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1340f32d5d43SBarry Smith         r1 = 0.0;
1341f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1342f32d5d43SBarry Smith         aj = a->j + a->i[i];
1343f32d5d43SBarry Smith         aa = a->a + a->i[i];
1344f32d5d43SBarry Smith 
1345f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1346f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1347f32d5d43SBarry Smith         }
1348a2ea699eSBarry Smith         c[colam + i] += r1;
1349f32d5d43SBarry Smith       }
1350f32d5d43SBarry Smith       b1 += bm;
1351f32d5d43SBarry Smith     }
13524324174eSBarry Smith   }
1353dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13548c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
13558c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13569123193aSHong Zhang   PetscFunctionReturn(0);
13579123193aSHong Zhang }
1358b1683b59SHong Zhang 
1359b1683b59SHong Zhang #undef __FUNCT__
1360b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1361b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1362c8db22f6SHong Zhang {
1363c8db22f6SHong Zhang   PetscErrorCode ierr;
13642f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13652f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13662f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13672f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13682f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13692f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1370c8db22f6SHong Zhang 
1371c8db22f6SHong Zhang   PetscFunctionBegin;
13722f5f1f90SHong Zhang   btval_den=btdense->v;
13732f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13742f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13752f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13762f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1377d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13782f5f1f90SHong Zhang       btcol = bj + bi[col];
13792f5f1f90SHong Zhang       btval = ba + bi[col];
13802f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1381d2853540SHong Zhang       for (j=0; j<anz; j++) {
13822f5f1f90SHong Zhang         brow            = btcol[j];
13832f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1384c8db22f6SHong Zhang       }
1385c8db22f6SHong Zhang     }
13862f5f1f90SHong Zhang     btval_den += m;
1387c8db22f6SHong Zhang   }
1388c8db22f6SHong Zhang   PetscFunctionReturn(0);
1389c8db22f6SHong Zhang }
1390c8db22f6SHong Zhang 
1391c8db22f6SHong Zhang #undef __FUNCT__
1392b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1393b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13948972f759SHong Zhang {
13958972f759SHong Zhang   PetscErrorCode ierr;
1396b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
13972f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1398b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
13992f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
14008972f759SHong Zhang 
14018972f759SHong Zhang   PetscFunctionBegin;
14020298fd71SBarry Smith   ierr   = MatGetLocalSize(Csp,&m,NULL);CHKERRQ(ierr);
1403a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1404b2d2b04fSHong Zhang   cp_den = ca_den;
14052f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
14062f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
14072f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
14082f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
14092f5f1f90SHong Zhang     for (l=0; l<nrows; l++) {
14102f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1411b2d2b04fSHong Zhang     }
1412b2d2b04fSHong Zhang     cp_den += m;
1413b2d2b04fSHong Zhang   }
1414a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
14158972f759SHong Zhang   PetscFunctionReturn(0);
14168972f759SHong Zhang }
14178972f759SHong Zhang 
1418e99be685SHong Zhang /*
1419e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1420e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1421e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1422e99be685SHong Zhang  */
1423e99be685SHong Zhang #undef __FUNCT__
1424e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
14251a83f524SJed 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)
1426e99be685SHong Zhang {
1427e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1428e99be685SHong Zhang   PetscErrorCode ierr;
1429e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1430e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1431e99be685SHong Zhang   PetscInt       *cspidx;
1432e99be685SHong Zhang 
1433e99be685SHong Zhang   PetscFunctionBegin;
1434e99be685SHong Zhang   *nn = n;
1435e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1436e99be685SHong Zhang   if (symmetric) {
1437ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
14381a83f524SJed Brown     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr);
1439e99be685SHong Zhang   } else {
1440e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1441e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1442e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1443e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1444e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1445e99be685SHong Zhang     jj   = a->j;
1446e99be685SHong Zhang     for (i=0; i<nz; i++) {
1447e99be685SHong Zhang       collengths[jj[i]]++;
1448e99be685SHong Zhang     }
1449e99be685SHong Zhang     cia[0] = oshift;
1450e99be685SHong Zhang     for (i=0; i<n; i++) {
1451e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1452e99be685SHong Zhang     }
1453e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1454e99be685SHong Zhang     jj   = a->j;
1455e99be685SHong Zhang     for (row=0; row<m; row++) {
1456e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1457e99be685SHong Zhang       for (i=0; i<mr; i++) {
1458e99be685SHong Zhang         col = *jj++;
14592205254eSKarl Rupp 
1460e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1461e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift]  = row + oshift;
1462e99be685SHong Zhang       }
1463e99be685SHong Zhang     }
1464e99be685SHong Zhang     ierr   = PetscFree(collengths);CHKERRQ(ierr);
1465e99be685SHong Zhang     *ia    = cia; *ja = cja;
1466e99be685SHong Zhang     *spidx = cspidx;
1467e99be685SHong Zhang   }
1468e99be685SHong Zhang   PetscFunctionReturn(0);
1469e99be685SHong Zhang }
1470e99be685SHong Zhang 
1471e99be685SHong Zhang #undef __FUNCT__
1472e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
14731a83f524SJed 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)
1474e99be685SHong Zhang {
1475e99be685SHong Zhang   PetscErrorCode ierr;
1476e99be685SHong Zhang 
1477e99be685SHong Zhang   PetscFunctionBegin;
1478e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1479e99be685SHong Zhang 
1480e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1481e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1482e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1483e99be685SHong Zhang   PetscFunctionReturn(0);
1484e99be685SHong Zhang }
1485e99be685SHong Zhang 
14868972f759SHong Zhang #undef __FUNCT__
1487b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1488b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1489b1683b59SHong Zhang {
1490b1683b59SHong Zhang   PetscErrorCode ierr;
14911a83f524SJed Brown   PetscInt       i,n,nrows,N,j,k,m,ncols,col,cm;
14921a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1493b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1494b1683b59SHong Zhang   IS             *isa;
1495b1683b59SHong Zhang   PetscBool      flg1,flg2;
1496b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1497e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1498d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1499b1683b59SHong Zhang 
1500b1683b59SHong Zhang   PetscFunctionBegin;
1501b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1502e99be685SHong Zhang 
1503b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1504251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1505251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1506b1683b59SHong Zhang   if (flg1 || flg2) {
1507b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1508b1683b59SHong Zhang   }
1509b1683b59SHong Zhang 
1510b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1511b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1512b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1513b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1514b1683b59SHong Zhang   c->rstart = 0;
1515b1683b59SHong Zhang 
1516b1683b59SHong Zhang   c->ncolors = nis;
1517b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1518b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1519d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1520d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
15212205254eSKarl Rupp 
1522d2853540SHong Zhang   colorforrow[0]    = 0;
1523d2853540SHong Zhang   rows_i            = rows;
1524d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1525b1683b59SHong Zhang 
1526d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1527d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
15282205254eSKarl Rupp 
1529d2853540SHong Zhang   colorforcol[0] = 0;
1530d2853540SHong Zhang   columns_i      = columns;
1531d2853540SHong Zhang 
15320298fd71SBarry Smith   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); /* column-wise storage of mat */
1533b1683b59SHong Zhang 
1534b28d80bdSHong Zhang   cm   = c->m;
1535b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1536e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1537b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1538b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1539b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
15402205254eSKarl Rupp 
1541b1683b59SHong Zhang     c->ncolumns[i] = n;
1542b1683b59SHong Zhang     if (n) {
1543d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1544b1683b59SHong Zhang     }
1545d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1546d2853540SHong Zhang     columns_i       += n;
1547b1683b59SHong Zhang 
1548b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1549e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1550e99be685SHong Zhang 
1551b1683b59SHong Zhang     /* loop over columns*/
1552b1683b59SHong Zhang     for (j=0; j<n; j++) {
1553b1683b59SHong Zhang       col     = is[j];
1554d2853540SHong Zhang       row_idx = cj + ci[col];
1555b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1556b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1557b1683b59SHong Zhang       for (k=0; k<m; k++) {
1558e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1559d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1560b1683b59SHong Zhang       }
1561b1683b59SHong Zhang     }
1562b1683b59SHong Zhang     /* count the number of hits */
1563b1683b59SHong Zhang     nrows = 0;
1564e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1565b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1566b1683b59SHong Zhang     }
1567b1683b59SHong Zhang     c->nrows[i]      = nrows;
1568d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1569d2853540SHong Zhang 
1570b1683b59SHong Zhang     nrows = 0;
1571e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1572b1683b59SHong Zhang       if (rowhit[j]) {
1573d2853540SHong Zhang         rows_i[nrows]            = j;
1574e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1575b1683b59SHong Zhang         nrows++;
1576b1683b59SHong Zhang       }
1577b1683b59SHong Zhang     }
1578b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1579d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1580b1683b59SHong Zhang   }
15810298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1582b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1583b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1584d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1585d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1586d2853540SHong Zhang #endif
1587b28d80bdSHong Zhang 
1588d2853540SHong Zhang   c->colorforrow     = colorforrow;
1589d2853540SHong Zhang   c->rows            = rows;
1590d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1591d2853540SHong Zhang   c->colorforcol     = colorforcol;
1592d2853540SHong Zhang   c->columns         = columns;
15932205254eSKarl Rupp 
1594f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1595b1683b59SHong Zhang   PetscFunctionReturn(0);
1596b1683b59SHong Zhang }
1597