xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 0ced3a2b8485496f9daf05a2cd5bf7bba7af8937)
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>
9*0ced3a2bSJed Brown #include <../src/mat/utils/petscheap.h>
10c6db04a5SJed Brown #include <petscbt.h>
11c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
12e9e4536cSHong Zhang /*
13e9e4536cSHong Zhang #define DEBUG_MATMATMULT
14e9e4536cSHong Zhang  */
15ec01deb9SMatthew Knepley EXTERN_C_BEGIN
166284ec50SHong Zhang #undef __FUNCT__
175c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1838baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1938baddfdSBarry Smith {
20dfbe8321SBarry Smith   PetscErrorCode ierr;
21*0ced3a2bSJed Brown   PetscBool      scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE;
225c66b693SKris Buschelman 
235c66b693SKris Buschelman   PetscFunctionBegin;
2426be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
25152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
26152983bfSHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
273c50cad2SHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,PETSC_NULL);CHKERRQ(ierr);
28*0ced3a2bSJed Brown     ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,PETSC_NULL);CHKERRQ(ierr);
29d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
30d8bbc50fSBarry Smith     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
313c50cad2SHong Zhang     if (scalable_fast){
323c50cad2SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
33aacf854cSBarry Smith     } else if (scalable){
34aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
35*0ced3a2bSJed Brown     } else if (heap) {
36*0ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
3725023028SHong Zhang     } else {
3826be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3925023028SHong Zhang     }
40598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
4126be0446SHong Zhang   }
42598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
4301e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
44598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
455c66b693SKris Buschelman   PetscFunctionReturn(0);
465c66b693SKris Buschelman }
47ec01deb9SMatthew Knepley EXTERN_C_END
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;
60fb908031SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_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);
77fb908031SHong Zhang   current_space = free_space;
78fb908031SHong Zhang 
79fb908031SHong Zhang   /* Determine ci and cj */
80fb908031SHong Zhang   for (i=0; i<am; i++) {
81fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
82fb908031SHong Zhang     aj   = a->j + ai[i];
83fb908031SHong Zhang     for (j=0; j<anzi; j++){
84fb908031SHong Zhang       brow = aj[j];
85fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
86fb908031SHong Zhang       bj   = b->j + bi[brow];
87fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
88fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
89fb908031SHong Zhang     }
90fb908031SHong Zhang     cnzi = lnk[0];
91fb908031SHong Zhang 
92fb908031SHong Zhang     /* If free space is not available, make more free space */
93fb908031SHong Zhang     /* Double the amount of total space in the list */
94fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
95fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
96fb908031SHong Zhang       ndouble++;
97fb908031SHong Zhang     }
98fb908031SHong Zhang 
99fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
100fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
101fb908031SHong Zhang     current_space->array           += cnzi;
102fb908031SHong Zhang     current_space->local_used      += cnzi;
103fb908031SHong Zhang     current_space->local_remaining -= cnzi;
104fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
105fb908031SHong Zhang   }
106fb908031SHong Zhang 
107fb908031SHong Zhang   /* Column indices are in the list of free space */
108fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
109fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
110fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
111fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
112fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
113b561aa9dSHong Zhang 
11426be0446SHong Zhang   /* put together the new symbolic matrix */
115aa1aec99SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr);
116a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
117a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
11858c24d83SHong Zhang 
11958c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
12058c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
12158c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
122aa1aec99SHong Zhang   c->free_a  = PETSC_FALSE;
123e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
12458c24d83SHong Zhang   c->nonew   = 0;
125002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1260b7e3e3dSHong Zhang 
1277212da7cSHong Zhang   /* set MatInfo */
1287212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
129f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1307212da7cSHong Zhang   c->maxnz                     = ci[am];
1317212da7cSHong Zhang   c->nz                        = ci[am];
132fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1337212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1347212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1357212da7cSHong Zhang 
1367212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1377212da7cSHong Zhang   if (ci[am]) {
138fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
139f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
140f2b054eeSHong Zhang   } else {
141f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
142be0fcf8dSHong Zhang   }
143f2b054eeSHong Zhang #endif
14458c24d83SHong Zhang   PetscFunctionReturn(0);
14558c24d83SHong Zhang }
146d50806bdSBarry Smith 
14726be0446SHong Zhang #undef __FUNCT__
14826be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
149dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
150d50806bdSBarry Smith {
151dfbe8321SBarry Smith   PetscErrorCode ierr;
152d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
153d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
154d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
155d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
15638baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
157c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
158519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
159aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
160aa1aec99SHong Zhang   PetscScalar    *ab_dense;
161d50806bdSBarry Smith 
162d50806bdSBarry Smith   PetscFunctionBegin;
163aa1aec99SHong Zhang   /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */
164aa1aec99SHong Zhang   if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
165aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
166aa1aec99SHong Zhang     c->a      = ca;
167aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
168aa1aec99SHong Zhang 
169aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
170aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
171aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
172aa1aec99SHong Zhang   } else {
173aa1aec99SHong Zhang     ca       = c->a;
174aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
175aa1aec99SHong Zhang   }
176aa1aec99SHong Zhang 
177c124e916SHong Zhang   /* clean old values in C */
178c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
179d50806bdSBarry Smith   /* Traverse A row-wise. */
180d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
181d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
182d50806bdSBarry Smith   for (i=0; i<am; i++) {
183d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
184d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
185519eb980SHong Zhang       brow = aj[j];
186d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
187d50806bdSBarry Smith       bjj  = bj + bi[brow];
188d50806bdSBarry Smith       baj  = ba + bi[brow];
189519eb980SHong Zhang       /* perform dense axpy */
19036ec6d2dSHong Zhang       valtmp = aa[j];
191519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
19236ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
193519eb980SHong Zhang       }
194519eb980SHong Zhang       flops += 2*bnzi;
195519eb980SHong Zhang     }
196c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
197c58ca2e3SHong Zhang 
198c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
199519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
200519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
201519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
202519eb980SHong Zhang     }
203519eb980SHong Zhang     flops += cnzi;
204519eb980SHong Zhang     cj += cnzi; ca += cnzi;
205519eb980SHong Zhang   }
206c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
207c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
208c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
209c58ca2e3SHong Zhang   PetscFunctionReturn(0);
210c58ca2e3SHong Zhang }
211c58ca2e3SHong Zhang 
212c58ca2e3SHong Zhang #undef __FUNCT__
21325023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
21425023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
215c58ca2e3SHong Zhang {
216c58ca2e3SHong Zhang   PetscErrorCode ierr;
217c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
218c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
219c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
220c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
221c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
222c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
223c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
22436ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
225c58ca2e3SHong Zhang   PetscInt       nextb;
226c58ca2e3SHong Zhang 
227c58ca2e3SHong Zhang   PetscFunctionBegin;
228e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
229971236abSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr);
230e9e4536cSHong Zhang #endif
231c58ca2e3SHong Zhang   /* clean old values in C */
232c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
233c58ca2e3SHong Zhang   /* Traverse A row-wise. */
234c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
235c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
236519eb980SHong Zhang   for (i=0;i<am;i++) {
237519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
238519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
239519eb980SHong Zhang     for (j=0;j<anzi;j++) {
240519eb980SHong Zhang       brow = aj[j];
241519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
242519eb980SHong Zhang       bjj  = bj + bi[brow];
243519eb980SHong Zhang       baj  = ba + bi[brow];
244519eb980SHong Zhang       /* perform sparse axpy */
24536ec6d2dSHong Zhang       valtmp = aa[j];
246c124e916SHong Zhang       nextb  = 0;
247c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
248c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
24936ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
250c124e916SHong Zhang         }
251d50806bdSBarry Smith       }
252d50806bdSBarry Smith       flops += 2*bnzi;
253d50806bdSBarry Smith     }
254519eb980SHong Zhang     aj += anzi; aa += anzi;
255519eb980SHong Zhang     cj += cnzi; ca += cnzi;
256d50806bdSBarry Smith   }
257c58ca2e3SHong Zhang 
258716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
261d50806bdSBarry Smith   PetscFunctionReturn(0);
262d50806bdSBarry Smith }
263bc011b1eSHong Zhang 
264e9e4536cSHong Zhang #undef __FUNCT__
2653c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2663c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
26725296bd5SBarry Smith {
26825296bd5SBarry Smith   PetscErrorCode     ierr;
26925296bd5SBarry Smith   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
27025296bd5SBarry Smith   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
2713c50cad2SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
27225296bd5SBarry Smith   MatScalar          *ca;
27325296bd5SBarry Smith   PetscReal          afill;
2743c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2753c50cad2SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
27625296bd5SBarry Smith 
27725296bd5SBarry Smith   PetscFunctionBegin;
2783c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2793c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2803c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2813c50cad2SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2823c50cad2SHong Zhang   ci[0] = 0;
2833c50cad2SHong Zhang 
2843c50cad2SHong Zhang   /* create and initialize a linked list */
2853c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2863c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2873c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
2883c50cad2SHong Zhang 
2893c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
2903c50cad2SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
2913c50cad2SHong Zhang   current_space = free_space;
2923c50cad2SHong Zhang 
2933c50cad2SHong Zhang   /* Determine ci and cj */
2943c50cad2SHong Zhang   for (i=0; i<am; i++) {
2953c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
2963c50cad2SHong Zhang     aj   = a->j + ai[i];
2973c50cad2SHong Zhang     for (j=0; j<anzi; j++){
2983c50cad2SHong Zhang       brow = aj[j];
2993c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3003c50cad2SHong Zhang       bj   = b->j + bi[brow];
3013c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3023c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3033c50cad2SHong Zhang     }
3043c50cad2SHong Zhang     cnzi = lnk[1];
3053c50cad2SHong Zhang 
3063c50cad2SHong Zhang     /* If free space is not available, make more free space */
3073c50cad2SHong Zhang     /* Double the amount of total space in the list */
3083c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3093c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3103c50cad2SHong Zhang       ndouble++;
3113c50cad2SHong Zhang     }
3123c50cad2SHong Zhang 
3133c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3143c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3153c50cad2SHong Zhang     current_space->array           += cnzi;
3163c50cad2SHong Zhang     current_space->local_used      += cnzi;
3173c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3183c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3193c50cad2SHong Zhang   }
3203c50cad2SHong Zhang 
3213c50cad2SHong Zhang   /* Column indices are in the list of free space */
3223c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3233c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3243c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3253c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3263c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
32725296bd5SBarry Smith 
32825296bd5SBarry Smith   /* Allocate space for ca */
32925296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
33025296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
33125296bd5SBarry Smith 
33225296bd5SBarry Smith   /* put together the new symbolic matrix */
33325296bd5SBarry Smith   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
334a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
335a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
33625296bd5SBarry Smith 
33725296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
33825296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
33925296bd5SBarry Smith   c = (Mat_SeqAIJ *)((*C)->data);
34025296bd5SBarry Smith   c->free_a   = PETSC_TRUE;
34125296bd5SBarry Smith   c->free_ij  = PETSC_TRUE;
34225296bd5SBarry Smith   c->nonew    = 0;
34325296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
34425296bd5SBarry Smith 
34525296bd5SBarry Smith   /* set MatInfo */
34625296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
34725296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
34825296bd5SBarry Smith   c->maxnz                     = ci[am];
34925296bd5SBarry Smith   c->nz                        = ci[am];
3503c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
35125296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
35225296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
35325296bd5SBarry Smith 
35425296bd5SBarry Smith #if defined(PETSC_USE_INFO)
35525296bd5SBarry Smith   if (ci[am]) {
3563c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
35725296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
35825296bd5SBarry Smith   } else {
35925296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
36025296bd5SBarry Smith   }
36125296bd5SBarry Smith #endif
36225296bd5SBarry Smith   PetscFunctionReturn(0);
36325296bd5SBarry Smith }
36425296bd5SBarry Smith 
36525296bd5SBarry Smith 
36625296bd5SBarry Smith #undef __FUNCT__
36725023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
36825023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
369e9e4536cSHong Zhang {
370e9e4536cSHong Zhang   PetscErrorCode     ierr;
371e9e4536cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
372bf9555e6SHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
37325c41797SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
374e9e4536cSHong Zhang   MatScalar          *ca;
375e9e4536cSHong Zhang   PetscReal          afill;
376bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
377bd958071SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
378e9e4536cSHong Zhang 
379e9e4536cSHong Zhang   PetscFunctionBegin;
380bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
381bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
382bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
383bd958071SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
384bd958071SHong Zhang   ci[0] = 0;
385bd958071SHong Zhang 
386bd958071SHong Zhang   /* create and initialize a linked list */
387bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
388bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
389bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
390bd958071SHong Zhang 
391bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
392bd958071SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
393bd958071SHong Zhang   current_space = free_space;
394bd958071SHong Zhang 
395bd958071SHong Zhang   /* Determine ci and cj */
396bd958071SHong Zhang   for (i=0; i<am; i++) {
397bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
398bd958071SHong Zhang     aj   = a->j + ai[i];
399bd958071SHong Zhang     for (j=0; j<anzi; j++){
400bd958071SHong Zhang       brow = aj[j];
401bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
402bd958071SHong Zhang       bj   = b->j + bi[brow];
403bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
404bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
405bd958071SHong Zhang     }
406bd958071SHong Zhang     cnzi = lnk[0];
407bd958071SHong Zhang 
408bd958071SHong Zhang     /* If free space is not available, make more free space */
409bd958071SHong Zhang     /* Double the amount of total space in the list */
410bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
411bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
412bd958071SHong Zhang       ndouble++;
413bd958071SHong Zhang     }
414bd958071SHong Zhang 
415bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
416bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
417bd958071SHong Zhang     current_space->array           += cnzi;
418bd958071SHong Zhang     current_space->local_used      += cnzi;
419bd958071SHong Zhang     current_space->local_remaining -= cnzi;
420bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
421bd958071SHong Zhang   }
422bd958071SHong Zhang 
423bd958071SHong Zhang   /* Column indices are in the list of free space */
424bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
425bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
426bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
427bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
428bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
429e9e4536cSHong Zhang 
430e9e4536cSHong Zhang   /* Allocate space for ca */
431bd958071SHong Zhang   /*-----------------------*/
432e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
433e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
434e9e4536cSHong Zhang 
435e9e4536cSHong Zhang   /* put together the new symbolic matrix */
436e9e4536cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
437a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
438a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
439e9e4536cSHong Zhang 
440e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
441e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
442e9e4536cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
443e9e4536cSHong Zhang   c->free_a   = PETSC_TRUE;
444e9e4536cSHong Zhang   c->free_ij  = PETSC_TRUE;
445e9e4536cSHong Zhang   c->nonew    = 0;
44625023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
447e9e4536cSHong Zhang 
448e9e4536cSHong Zhang   /* set MatInfo */
449e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
450e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
451e9e4536cSHong Zhang   c->maxnz                     = ci[am];
452e9e4536cSHong Zhang   c->nz                        = ci[am];
453bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
454e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
455e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
456e9e4536cSHong Zhang 
457e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
458e9e4536cSHong Zhang   if (ci[am]) {
459bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
460e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
461e9e4536cSHong Zhang   } else {
462e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
463e9e4536cSHong Zhang   }
464e9e4536cSHong Zhang #endif
465e9e4536cSHong Zhang   PetscFunctionReturn(0);
466e9e4536cSHong Zhang }
467e9e4536cSHong Zhang 
468*0ced3a2bSJed Brown #undef __FUNCT__
469*0ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
470*0ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
471*0ced3a2bSJed Brown {
472*0ced3a2bSJed Brown   PetscErrorCode     ierr;
473*0ced3a2bSJed Brown   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
474*0ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
475*0ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
476*0ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
477*0ced3a2bSJed Brown   MatScalar          *ca;
478*0ced3a2bSJed Brown   PetscReal          afill;
479*0ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
480*0ced3a2bSJed Brown   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
481*0ced3a2bSJed Brown   PetscHeap          h;
482*0ced3a2bSJed Brown 
483*0ced3a2bSJed Brown   PetscFunctionBegin;
484*0ced3a2bSJed Brown   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
485*0ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
486*0ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
487*0ced3a2bSJed Brown   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
488*0ced3a2bSJed Brown   ci[0] = 0;
489*0ced3a2bSJed Brown 
490*0ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
491*0ced3a2bSJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
492*0ced3a2bSJed Brown   current_space = free_space;
493*0ced3a2bSJed Brown 
494*0ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
495*0ced3a2bSJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
496*0ced3a2bSJed Brown 
497*0ced3a2bSJed Brown   /* Determine ci and cj */
498*0ced3a2bSJed Brown   for (i=0; i<am; i++) {
499*0ced3a2bSJed 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 */
500*0ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
501*0ced3a2bSJed Brown     ci[i+1] = ci[i];
502*0ced3a2bSJed Brown     /* Populate the min heap */
503*0ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
504*0ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
505*0ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
506*0ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
507*0ced3a2bSJed Brown       }
508*0ced3a2bSJed Brown     }
509*0ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
510*0ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
511*0ced3a2bSJed Brown     while (j >= 0) {
512*0ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
513*0ced3a2bSJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
514*0ced3a2bSJed Brown         ndouble++;
515*0ced3a2bSJed Brown       }
516*0ced3a2bSJed Brown       *(current_space->array++) = col;
517*0ced3a2bSJed Brown       current_space->local_used++;
518*0ced3a2bSJed Brown       current_space->local_remaining--;
519*0ced3a2bSJed Brown       ci[i+1]++;
520*0ced3a2bSJed Brown 
521*0ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
522*0ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
523*0ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
524*0ced3a2bSJed Brown         PetscInt j2,col2;
525*0ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
526*0ced3a2bSJed Brown         if (col2 != col) break;
527*0ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
528*0ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
529*0ced3a2bSJed Brown       }
530*0ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
531*0ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
532*0ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
533*0ced3a2bSJed Brown     }
534*0ced3a2bSJed Brown   }
535*0ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
536*0ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
537*0ced3a2bSJed Brown 
538*0ced3a2bSJed Brown   /* Column indices are in the list of free space */
539*0ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
540*0ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
541*0ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
542*0ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
543*0ced3a2bSJed Brown 
544*0ced3a2bSJed Brown   /* Allocate space for ca */
545*0ced3a2bSJed Brown   /*-----------------------*/
546*0ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(MatScalar),&ca);CHKERRQ(ierr);
547*0ced3a2bSJed Brown   ierr = PetscMemzero(ca,ci[am]*sizeof(MatScalar));CHKERRQ(ierr);
548*0ced3a2bSJed Brown 
549*0ced3a2bSJed Brown   /* put together the new symbolic matrix */
550*0ced3a2bSJed Brown   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
551*0ced3a2bSJed Brown   (*C)->rmap->bs = A->rmap->bs;
552*0ced3a2bSJed Brown   (*C)->cmap->bs = B->cmap->bs;
553*0ced3a2bSJed Brown 
554*0ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
555*0ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
556*0ced3a2bSJed Brown   c = (Mat_SeqAIJ *)((*C)->data);
557*0ced3a2bSJed Brown   c->free_a   = PETSC_TRUE;
558*0ced3a2bSJed Brown   c->free_ij  = PETSC_TRUE;
559*0ced3a2bSJed Brown   c->nonew    = 0;
560*0ced3a2bSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
561*0ced3a2bSJed Brown 
562*0ced3a2bSJed Brown   /* set MatInfo */
563*0ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
564*0ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
565*0ced3a2bSJed Brown   c->maxnz                     = ci[am];
566*0ced3a2bSJed Brown   c->nz                        = ci[am];
567*0ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
568*0ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
569*0ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
570*0ced3a2bSJed Brown 
571*0ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
572*0ced3a2bSJed Brown   if (ci[am]) {
573*0ced3a2bSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
574*0ced3a2bSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
575*0ced3a2bSJed Brown   } else {
576*0ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
577*0ced3a2bSJed Brown   }
578*0ced3a2bSJed Brown #endif
579*0ced3a2bSJed Brown   PetscFunctionReturn(0);
580*0ced3a2bSJed Brown }
581e9e4536cSHong Zhang 
582d2853540SHong Zhang /* This routine is not used. Should be removed! */
583bc011b1eSHong Zhang #undef __FUNCT__
5846fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
5856fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
5865df89d91SHong Zhang {
587bc011b1eSHong Zhang   PetscErrorCode ierr;
588bc011b1eSHong Zhang 
589bc011b1eSHong Zhang   PetscFunctionBegin;
590bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5916fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
592bc011b1eSHong Zhang   }
5936fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
594bc011b1eSHong Zhang   PetscFunctionReturn(0);
595bc011b1eSHong Zhang }
596bc011b1eSHong Zhang 
597bc011b1eSHong Zhang #undef __FUNCT__
598e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
599e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
6002128a86cSHong Zhang {
6012128a86cSHong Zhang   PetscErrorCode      ierr;
602e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
6032128a86cSHong Zhang 
6042128a86cSHong Zhang   PetscFunctionBegin;
605b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
6062128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
6072128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
6082128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
6092128a86cSHong Zhang   PetscFunctionReturn(0);
6102128a86cSHong Zhang }
6112128a86cSHong Zhang 
6122128a86cSHong Zhang #undef __FUNCT__
6132128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
6142128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
6152128a86cSHong Zhang {
6162128a86cSHong Zhang   PetscErrorCode      ierr;
6172128a86cSHong Zhang   PetscContainer      container;
618e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
6192128a86cSHong Zhang 
6202128a86cSHong Zhang   PetscFunctionBegin;
621e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
6222128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
6232128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
6242128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
6252128a86cSHong Zhang   if (A->ops->destroy) {
6262128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
6272128a86cSHong Zhang   }
628e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
6292128a86cSHong Zhang   PetscFunctionReturn(0);
6302128a86cSHong Zhang }
6312128a86cSHong Zhang 
6322128a86cSHong Zhang #undef __FUNCT__
6336fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
6346fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
635bc011b1eSHong Zhang {
6365df89d91SHong Zhang   PetscErrorCode      ierr;
63781d82fe4SHong Zhang   Mat                 Bt;
63881d82fe4SHong Zhang   PetscInt            *bti,*btj;
639e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
6402128a86cSHong Zhang   PetscContainer      container;
641d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
642d2853540SHong Zhang 
64381d82fe4SHong Zhang   PetscFunctionBegin;
644d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
64581d82fe4SHong Zhang    /* create symbolic Bt */
64681d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
64781d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
648c3e5699bSSatish Balay   Bt->rmap->bs = A->cmap->bs;
649c3e5699bSSatish Balay   Bt->cmap->bs = B->cmap->bs;
65081d82fe4SHong Zhang 
65181d82fe4SHong Zhang   /* get symbolic C=A*Bt */
65281d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
65381d82fe4SHong Zhang 
6542128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
655e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
6562128a86cSHong Zhang 
6572128a86cSHong Zhang   /* attach the supporting struct to C */
6582128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
6592128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
660e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
661e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
6622128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
6632128a86cSHong Zhang 
6642128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
6652128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
6662128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
6672128a86cSHong Zhang 
668d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
669d2853540SHong Zhang   etime2 += tf - t0;
670d2853540SHong Zhang 
671f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
6722128a86cSHong Zhang   if (multtrans->usecoloring){
673b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
674b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
6752128a86cSHong Zhang     ISColoring           iscoloring;
6762128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
677d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
678e8354b3eSHong Zhang 
679e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
680502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
681d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
682d2853540SHong Zhang     etime0 += tf - t0;
683d2853540SHong Zhang 
684d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
685b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
6862128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
6872128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
688e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
689d2853540SHong Zhang     etime01 += tf - t0;
6902128a86cSHong Zhang 
691e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
6922128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
6932128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
6942128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
6952128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
6962128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
6972128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
6982128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
6992128a86cSHong Zhang 
7002128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
7012128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
7022128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
7032128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
7042128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
7052128a86cSHong Zhang     multtrans->ABt_den = C_dense;
706e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
707e8354b3eSHong Zhang     etime1 += tf - t0;
708f94ccd6cSHong Zhang 
709f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
7101ce71dffSSatish Balay     {
711f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
712f94ccd6cSHong Zhang     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));
713f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
7141ce71dffSSatish Balay     }
715f94ccd6cSHong Zhang #endif
7162128a86cSHong Zhang   }
717e99be685SHong Zhang   /* clean up */
718e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
719e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
7202128a86cSHong Zhang 
721f94ccd6cSHong Zhang 
722f94ccd6cSHong Zhang 
72381d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
72481d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
7251fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
7261fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
7271fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
7281fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
7291fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
7301fa1209cSHong Zhang   MatScalar          *ca;
7311fa1209cSHong Zhang   PetscBT            lnkbt;
7321fa1209cSHong Zhang   PetscReal          afill;
7335df89d91SHong Zhang 
7341fa1209cSHong Zhang   /* Allocate row pointer array ci  */
7351fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
7361fa1209cSHong Zhang   ci[0] = 0;
7371fa1209cSHong Zhang 
7381fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
7391fa1209cSHong Zhang   nlnk = bm+1;
7401fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
7411fa1209cSHong Zhang 
7421fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
7431fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
7441fa1209cSHong Zhang   current_space = free_space;
7451fa1209cSHong Zhang 
7461fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
7471fa1209cSHong Zhang   for (i=0; i<am; i++) {
7481fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
7491fa1209cSHong Zhang     cnzi = 0;
7501fa1209cSHong Zhang     acol = aj + ai[i];
7511fa1209cSHong Zhang     for (j=0; j<bm; j++){
7521fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
7531fa1209cSHong Zhang       bcol= bj + bi[j];
75481d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
7551fa1209cSHong Zhang       ka = 0; kb = 0;
7561fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
75781d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
75881d82fe4SHong Zhang         if (ka == anzi) break;
75981d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
76081d82fe4SHong Zhang         if (kb == bnzj) break;
7611fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
7621fa1209cSHong Zhang           index[0] = j;
7631fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
7641fa1209cSHong Zhang           cnzi++;
7651fa1209cSHong Zhang           break;
7661fa1209cSHong Zhang         }
7671fa1209cSHong Zhang       }
7681fa1209cSHong Zhang     }
7691fa1209cSHong Zhang 
7701fa1209cSHong Zhang     /* If free space is not available, make more free space */
7711fa1209cSHong Zhang     /* Double the amount of total space in the list */
7721fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
7731fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
7741fa1209cSHong Zhang       nspacedouble++;
7751fa1209cSHong Zhang     }
7761fa1209cSHong Zhang 
7771fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
7781fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
7791fa1209cSHong Zhang     current_space->array           += cnzi;
7801fa1209cSHong Zhang     current_space->local_used      += cnzi;
7811fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
7821fa1209cSHong Zhang 
7831fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
7841fa1209cSHong Zhang   }
7851fa1209cSHong Zhang 
7861fa1209cSHong Zhang 
7871fa1209cSHong Zhang   /* Column indices are in the list of free space.
7881fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
7891fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
7901fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
7911fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7921fa1209cSHong Zhang 
7931fa1209cSHong Zhang   /* Allocate space for ca */
7941fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
7951fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
7961fa1209cSHong Zhang 
7971fa1209cSHong Zhang   /* put together the new symbolic matrix */
7981fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
799a2f3521dSMark F. Adams   (*C)->rmap->bs = A->cmap->bs;
800a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
8011fa1209cSHong Zhang 
8021fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
8031fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
8041fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
8051fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
8061fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
8071fa1209cSHong Zhang   c->nonew    = 0;
8081fa1209cSHong Zhang 
8091fa1209cSHong Zhang   /* set MatInfo */
8101fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
8111fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
8121fa1209cSHong Zhang   c->maxnz                     = ci[am];
8131fa1209cSHong Zhang   c->nz                        = ci[am];
8141fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
8151fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
8161fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
8171fa1209cSHong Zhang 
8181fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
8191fa1209cSHong Zhang   if (ci[am]) {
8201fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
8216fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
8221fa1209cSHong Zhang   } else {
8231fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
8241fa1209cSHong Zhang   }
8251fa1209cSHong Zhang #endif
82681d82fe4SHong Zhang #endif
8275df89d91SHong Zhang   PetscFunctionReturn(0);
8285df89d91SHong Zhang }
8295df89d91SHong Zhang 
8306973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
8315df89d91SHong Zhang #undef __FUNCT__
8326fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
8336fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
8345df89d91SHong Zhang {
8355df89d91SHong Zhang   PetscErrorCode ierr;
8365df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
837e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
83881d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
8395df89d91SHong Zhang   PetscLogDouble flops=0.0;
840aa1aec99SHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval;
841e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
8422128a86cSHong Zhang   PetscContainer      container;
8436973769bSHong Zhang #if defined(USE_ARRAY)
8446973769bSHong Zhang   MatScalar      *spdot;
8456973769bSHong Zhang #endif
8465df89d91SHong Zhang 
8475df89d91SHong Zhang   PetscFunctionBegin;
848e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
8492128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
8502128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
8512128a86cSHong Zhang   if (multtrans->usecoloring){
852b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
853c8db22f6SHong Zhang     Mat                   Bt_dense;
854c8db22f6SHong Zhang     PetscInt              m,n;
855b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
856a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
857a0b95ffbSSatish Balay 
8582128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
859c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
860c8db22f6SHong Zhang 
861b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
8622128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
863b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
8642128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
865b2d2b04fSHong Zhang     etime0 += tf - t0;
866c8db22f6SHong Zhang 
867c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
8682128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
8692128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
8702128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
8712128a86cSHong Zhang     etime2 += tf - t0;
872c8db22f6SHong Zhang 
8732128a86cSHong Zhang     /* Recover C from C_dense */
8742128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
875b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
8762128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
8772128a86cSHong Zhang     etime1 += tf - t0;
878f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
879f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
880f94ccd6cSHong Zhang #endif
881c8db22f6SHong Zhang     PetscFunctionReturn(0);
882c8db22f6SHong Zhang   }
883c8db22f6SHong Zhang 
8846973769bSHong Zhang #if defined(USE_ARRAY)
8856973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
886e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
887e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
8886973769bSHong Zhang #endif
8896973769bSHong Zhang 
89081d82fe4SHong Zhang   /* clear old values in C */
891aa1aec99SHong Zhang   if (!c->a){
892aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
893aa1aec99SHong Zhang     c->a      = ca;
894aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
895aa1aec99SHong Zhang   } else {
896aa1aec99SHong Zhang     ca =  c->a;
897aa1aec99SHong Zhang   }
89881d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
8995df89d91SHong Zhang 
9001fa1209cSHong Zhang   for (i=0; i<cm; i++) {
90181d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
90281d82fe4SHong Zhang     acol = aj + ai[i];
9036973769bSHong Zhang     aval = aa + ai[i];
9041fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9051fa1209cSHong Zhang     ccol = cj + ci[i];
9066973769bSHong Zhang     cval = ca + ci[i];
9071fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
90881d82fe4SHong Zhang       brow = ccol[j];
90981d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
91081d82fe4SHong Zhang       bcol = bj + bi[brow];
9116973769bSHong Zhang       bval = ba + bi[brow];
9126973769bSHong Zhang 
91381d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
9146973769bSHong Zhang #if defined(USE_ARRAY)
9156973769bSHong Zhang       /* put ba to spdot array */
9166973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
9176973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
9186973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
9196973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
9206973769bSHong Zhang       }
9216973769bSHong Zhang       /* zero spdot array */
9226973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
9236973769bSHong Zhang #else
92481d82fe4SHong Zhang       nexta = 0; nextb = 0;
92581d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
92681d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
92781d82fe4SHong Zhang         if (nexta == anzi) break;
92881d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
92981d82fe4SHong Zhang         if (nextb == bnzj) break;
93081d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
9316973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
93281d82fe4SHong Zhang           nexta++; nextb++;
93381d82fe4SHong Zhang           flops += 2;
9341fa1209cSHong Zhang         }
9351fa1209cSHong Zhang       }
9366973769bSHong Zhang #endif
93781d82fe4SHong Zhang     }
93881d82fe4SHong Zhang   }
93981d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94081d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94181d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
9426973769bSHong Zhang #if defined(USE_ARRAY)
9436973769bSHong Zhang   ierr = PetscFree(spdot);
9446973769bSHong Zhang #endif
9455df89d91SHong Zhang   PetscFunctionReturn(0);
9465df89d91SHong Zhang }
9475df89d91SHong Zhang 
9485df89d91SHong Zhang #undef __FUNCT__
94975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
95075648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
9515df89d91SHong Zhang   PetscErrorCode ierr;
9525df89d91SHong Zhang 
9535df89d91SHong Zhang   PetscFunctionBegin;
9545df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
95575648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
9565df89d91SHong Zhang   }
95775648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
9585df89d91SHong Zhang   PetscFunctionReturn(0);
9595df89d91SHong Zhang }
9605df89d91SHong Zhang 
9615df89d91SHong Zhang #undef __FUNCT__
96275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
96375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
9645df89d91SHong Zhang {
965bc011b1eSHong Zhang   PetscErrorCode ierr;
966bc011b1eSHong Zhang   Mat            At;
96738baddfdSBarry Smith   PetscInt       *ati,*atj;
968bc011b1eSHong Zhang 
969bc011b1eSHong Zhang   PetscFunctionBegin;
970bc011b1eSHong Zhang   /* create symbolic At */
971bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
972d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
973a2f3521dSMark F. Adams   At->rmap->bs = A->cmap->bs;
974a2f3521dSMark F. Adams   At->cmap->bs = B->cmap->bs;
975bc011b1eSHong Zhang 
976bc011b1eSHong Zhang   /* get symbolic C=At*B */
977bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
978bc011b1eSHong Zhang 
979bc011b1eSHong Zhang   /* clean up */
9806bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
981bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
982bc011b1eSHong Zhang   PetscFunctionReturn(0);
983bc011b1eSHong Zhang }
984bc011b1eSHong Zhang 
985bc011b1eSHong Zhang #undef __FUNCT__
98675648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
98775648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
988bc011b1eSHong Zhang {
989bc011b1eSHong Zhang   PetscErrorCode ierr;
9900fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
991d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
992d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
993d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
994aa1aec99SHong Zhang   MatScalar      *aa=a->a,*ba,*ca,*caj;
995bc011b1eSHong Zhang 
996bc011b1eSHong Zhang   PetscFunctionBegin;
997aa1aec99SHong Zhang   if (!c->a){
998aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
999aa1aec99SHong Zhang     c->a      = ca;
1000aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1001aa1aec99SHong Zhang   } else {
1002aa1aec99SHong Zhang     ca = c->a;
1003aa1aec99SHong Zhang   }
1004bc011b1eSHong Zhang   /* clear old values in C */
1005bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1006bc011b1eSHong Zhang 
1007bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1008bc011b1eSHong Zhang   for (i=0;i<am;i++) {
1009bc011b1eSHong Zhang     bj   = b->j + bi[i];
1010bc011b1eSHong Zhang     ba   = b->a + bi[i];
1011bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1012bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1013bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1014bc011b1eSHong Zhang       nextb = 0;
10150fbc74f4SHong Zhang       crow  = *aj++;
1016bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1017bc011b1eSHong Zhang       caj   = ca + ci[crow];
1018bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1019bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
10200fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
10210fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1022bc011b1eSHong Zhang           nextb++;
1023bc011b1eSHong Zhang         }
1024bc011b1eSHong Zhang       }
1025bc011b1eSHong Zhang       flops += 2*bnzi;
10260fbc74f4SHong Zhang       aa++;
1027bc011b1eSHong Zhang     }
1028bc011b1eSHong Zhang   }
1029bc011b1eSHong Zhang 
1030bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1031bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1032bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1033bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1034bc011b1eSHong Zhang   PetscFunctionReturn(0);
1035bc011b1eSHong Zhang }
10369123193aSHong Zhang 
1037ec01deb9SMatthew Knepley EXTERN_C_BEGIN
10389123193aSHong Zhang #undef __FUNCT__
10399123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
10409123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10419123193aSHong Zhang {
10429123193aSHong Zhang   PetscErrorCode ierr;
10439123193aSHong Zhang 
10449123193aSHong Zhang   PetscFunctionBegin;
10459123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
10469123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
10479123193aSHong Zhang   }
10489123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
10499123193aSHong Zhang   PetscFunctionReturn(0);
10509123193aSHong Zhang }
1051ec01deb9SMatthew Knepley EXTERN_C_END
1052edd81eecSMatthew Knepley 
10539123193aSHong Zhang #undef __FUNCT__
10549123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
10559123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
10569123193aSHong Zhang {
10579123193aSHong Zhang   PetscErrorCode ierr;
10589123193aSHong Zhang 
10599123193aSHong Zhang   PetscFunctionBegin;
10605a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
10618cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
10629123193aSHong Zhang   PetscFunctionReturn(0);
10639123193aSHong Zhang }
10649123193aSHong Zhang 
10659123193aSHong Zhang #undef __FUNCT__
10669123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
10679123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
10689123193aSHong Zhang {
1069f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
10709123193aSHong Zhang   PetscErrorCode ierr;
1071dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1072dd6ea824SBarry Smith   MatScalar      *aa;
1073d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
1074f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
10759123193aSHong Zhang 
10769123193aSHong Zhang   PetscFunctionBegin;
1077f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1078e32f2f54SBarry 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);
1079e32f2f54SBarry 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);
1080e32f2f54SBarry 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);
1081f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
1082f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
1083f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1084f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
1085f32d5d43SBarry Smith     colam = col*am;
1086f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1087f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1088f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
1089f32d5d43SBarry Smith       aj  = a->j + a->i[i];
1090f32d5d43SBarry Smith       aa  = a->a + a->i[i];
1091f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1092f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
1093f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
1094f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
1095f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
10969123193aSHong Zhang       }
1097f32d5d43SBarry Smith       c[colam + i]       = r1;
1098f32d5d43SBarry Smith       c[colam + am + i]  = r2;
1099f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
1100f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
1101f32d5d43SBarry Smith     }
1102f32d5d43SBarry Smith     b1 += bm4;
1103f32d5d43SBarry Smith     b2 += bm4;
1104f32d5d43SBarry Smith     b3 += bm4;
1105f32d5d43SBarry Smith     b4 += bm4;
1106f32d5d43SBarry Smith   }
1107f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
1108f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1109f32d5d43SBarry Smith       r1 = 0.0;
1110f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
1111f32d5d43SBarry Smith       aj  = a->j + a->i[i];
1112f32d5d43SBarry Smith       aa  = a->a + a->i[i];
1113f32d5d43SBarry Smith 
1114f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1115f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
1116f32d5d43SBarry Smith       }
1117f32d5d43SBarry Smith       c[col*am + i]     = r1;
1118f32d5d43SBarry Smith     }
1119f32d5d43SBarry Smith     b1 += bm;
1120f32d5d43SBarry Smith   }
1121dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
1122f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1123f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
1124f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1125f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1126f32d5d43SBarry Smith   PetscFunctionReturn(0);
1127f32d5d43SBarry Smith }
1128f32d5d43SBarry Smith 
1129f32d5d43SBarry Smith /*
11304324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1131f32d5d43SBarry Smith */
1132f32d5d43SBarry Smith #undef __FUNCT__
1133f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1134f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1135f32d5d43SBarry Smith {
1136f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1137f32d5d43SBarry Smith   PetscErrorCode ierr;
1138dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1139dd6ea824SBarry Smith   MatScalar      *aa;
1140d0f46423SBarry 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;
11414324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1142f32d5d43SBarry Smith 
1143f32d5d43SBarry Smith   PetscFunctionBegin;
1144f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1145f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
1146f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
1147f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
11484324174eSBarry Smith 
11494324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
11504324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
11514324174eSBarry Smith       colam = col*am;
11524324174eSBarry Smith       arm   = a->compressedrow.nrows;
11534324174eSBarry Smith       ii    = a->compressedrow.i;
11544324174eSBarry Smith       ridx  = a->compressedrow.rindex;
11554324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
11564324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
11574324174eSBarry Smith 	n   = ii[i+1] - ii[i];
11584324174eSBarry Smith 	aj  = a->j + ii[i];
11594324174eSBarry Smith 	aa  = a->a + ii[i];
11604324174eSBarry Smith 	for (j=0; j<n; j++) {
11614324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
11624324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
11634324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
11644324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
11654324174eSBarry Smith 	}
11664324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
11674324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
11684324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
11694324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
11704324174eSBarry Smith       }
11714324174eSBarry Smith       b1 += bm4;
11724324174eSBarry Smith       b2 += bm4;
11734324174eSBarry Smith       b3 += bm4;
11744324174eSBarry Smith       b4 += bm4;
11754324174eSBarry Smith     }
11764324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
11774324174eSBarry Smith       colam = col*am;
11784324174eSBarry Smith       arm   = a->compressedrow.nrows;
11794324174eSBarry Smith       ii    = a->compressedrow.i;
11804324174eSBarry Smith       ridx  = a->compressedrow.rindex;
11814324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
11824324174eSBarry Smith 	r1 = 0.0;
11834324174eSBarry Smith 	n   = ii[i+1] - ii[i];
11844324174eSBarry Smith 	aj  = a->j + ii[i];
11854324174eSBarry Smith 	aa  = a->a + ii[i];
11864324174eSBarry Smith 
11874324174eSBarry Smith 	for (j=0; j<n; j++) {
11884324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
11894324174eSBarry Smith 	}
11904324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
11914324174eSBarry Smith       }
11924324174eSBarry Smith       b1 += bm;
11934324174eSBarry Smith     }
11944324174eSBarry Smith   } else {
1195f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
1196f32d5d43SBarry Smith       colam = col*am;
1197f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1198f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
1199f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1200f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1201f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1202f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1203f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
1204f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
1205f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
1206f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
1207f32d5d43SBarry Smith 	}
1208f32d5d43SBarry Smith 	c[colam + i]       += r1;
1209f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
1210f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
1211f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
1212f32d5d43SBarry Smith       }
1213f32d5d43SBarry Smith       b1 += bm4;
1214f32d5d43SBarry Smith       b2 += bm4;
1215f32d5d43SBarry Smith       b3 += bm4;
1216f32d5d43SBarry Smith       b4 += bm4;
1217f32d5d43SBarry Smith     }
1218f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
1219f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1220f32d5d43SBarry Smith 	r1 = 0.0;
1221f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1222f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1223f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1224f32d5d43SBarry Smith 
1225f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1226f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
1227f32d5d43SBarry Smith 	}
1228f32d5d43SBarry Smith 	c[col*am + i]     += r1;
1229f32d5d43SBarry Smith       }
1230f32d5d43SBarry Smith       b1 += bm;
1231f32d5d43SBarry Smith     }
12324324174eSBarry Smith   }
1233dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
1234f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1235f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
12369123193aSHong Zhang   PetscFunctionReturn(0);
12379123193aSHong Zhang }
1238b1683b59SHong Zhang 
1239b1683b59SHong Zhang #undef __FUNCT__
1240b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1241b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1242c8db22f6SHong Zhang {
1243c8db22f6SHong Zhang   PetscErrorCode ierr;
12442f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
12452f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
12462f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
12472f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
12482f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
12492f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1250c8db22f6SHong Zhang 
1251c8db22f6SHong Zhang   PetscFunctionBegin;
12522f5f1f90SHong Zhang   btval_den=btdense->v;
12532f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
12542f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
12552f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
12562f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1257d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
12582f5f1f90SHong Zhang       btcol = bj + bi[col];
12592f5f1f90SHong Zhang       btval = ba + bi[col];
12602f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1261d2853540SHong Zhang       for (j=0; j<anz; j++){
12622f5f1f90SHong Zhang         brow            = btcol[j];
12632f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1264c8db22f6SHong Zhang       }
1265c8db22f6SHong Zhang     }
12662f5f1f90SHong Zhang     btval_den += m;
1267c8db22f6SHong Zhang   }
1268c8db22f6SHong Zhang   PetscFunctionReturn(0);
1269c8db22f6SHong Zhang }
1270c8db22f6SHong Zhang 
1271c8db22f6SHong Zhang #undef __FUNCT__
1272b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1273b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
12748972f759SHong Zhang {
12758972f759SHong Zhang   PetscErrorCode ierr;
1276b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
12772f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1278b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
12792f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
12808972f759SHong Zhang 
12818972f759SHong Zhang   PetscFunctionBegin;
12828972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
1283b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
1284b2d2b04fSHong Zhang   cp_den = ca_den;
12852f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
12862f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
12872f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
12882f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
12892f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
12902f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1291b2d2b04fSHong Zhang     }
1292b2d2b04fSHong Zhang     cp_den += m;
1293b2d2b04fSHong Zhang   }
1294b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
12958972f759SHong Zhang   PetscFunctionReturn(0);
12968972f759SHong Zhang }
12978972f759SHong Zhang 
1298e99be685SHong Zhang /*
1299e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1300e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1301e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1302e99be685SHong Zhang  */
1303e99be685SHong Zhang #undef __FUNCT__
1304e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
1305e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1306e99be685SHong Zhang {
1307e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1308e99be685SHong Zhang   PetscErrorCode ierr;
1309e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1310e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1311e99be685SHong Zhang   PetscInt       *cspidx;
1312e99be685SHong Zhang 
1313e99be685SHong Zhang   PetscFunctionBegin;
1314e99be685SHong Zhang   *nn = n;
1315e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1316e99be685SHong Zhang   if (symmetric) {
1317e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
1318e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
1319e99be685SHong Zhang   } else {
1320e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1321e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1322e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1323e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1324e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1325e99be685SHong Zhang     jj = a->j;
1326e99be685SHong Zhang     for (i=0; i<nz; i++) {
1327e99be685SHong Zhang       collengths[jj[i]]++;
1328e99be685SHong Zhang     }
1329e99be685SHong Zhang     cia[0] = oshift;
1330e99be685SHong Zhang     for (i=0; i<n; i++) {
1331e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1332e99be685SHong Zhang     }
1333e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1334e99be685SHong Zhang     jj   = a->j;
1335e99be685SHong Zhang     for (row=0; row<m; row++) {
1336e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1337e99be685SHong Zhang       for (i=0; i<mr; i++) {
1338e99be685SHong Zhang         col = *jj++;
1339e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1340e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1341e99be685SHong Zhang       }
1342e99be685SHong Zhang     }
1343e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1344e99be685SHong Zhang     *ia = cia; *ja = cja;
1345e99be685SHong Zhang     *spidx = cspidx;
1346e99be685SHong Zhang   }
1347e99be685SHong Zhang   PetscFunctionReturn(0);
1348e99be685SHong Zhang }
1349e99be685SHong Zhang 
1350e99be685SHong Zhang #undef __FUNCT__
1351e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1352e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1353e99be685SHong Zhang {
1354e99be685SHong Zhang   PetscErrorCode ierr;
1355e99be685SHong Zhang 
1356e99be685SHong Zhang   PetscFunctionBegin;
1357e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1358e99be685SHong Zhang 
1359e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1360e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1361e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1362e99be685SHong Zhang   PetscFunctionReturn(0);
1363e99be685SHong Zhang }
1364e99be685SHong Zhang 
13658972f759SHong Zhang #undef __FUNCT__
1366b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1367b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1368b1683b59SHong Zhang {
1369b1683b59SHong Zhang   PetscErrorCode ierr;
1370d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1371b1683b59SHong Zhang   const PetscInt *is;
1372b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1373b1683b59SHong Zhang   IS             *isa;
1374b1683b59SHong Zhang   PetscBool      done;
1375b1683b59SHong Zhang   PetscBool      flg1,flg2;
1376b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1377e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1378d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1379b1683b59SHong Zhang 
1380b1683b59SHong Zhang   PetscFunctionBegin;
1381b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1382e99be685SHong Zhang 
1383b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1384251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1385251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1386b1683b59SHong Zhang   if (flg1 || flg2) {
1387b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1388b1683b59SHong Zhang   }
1389b1683b59SHong Zhang 
1390b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1391b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1392b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1393b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1394b1683b59SHong Zhang   c->rstart = 0;
1395b1683b59SHong Zhang 
1396b1683b59SHong Zhang   c->ncolors = nis;
1397b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1398b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1399d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1400d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1401d2853540SHong Zhang   colorforrow[0]    = 0;
1402d2853540SHong Zhang   rows_i            = rows;
1403d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1404b1683b59SHong Zhang 
1405d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1406d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1407d2853540SHong Zhang   colorforcol[0] = 0;
1408d2853540SHong Zhang   columns_i      = columns;
1409d2853540SHong Zhang 
1410e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1411b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1412b1683b59SHong Zhang 
1413b28d80bdSHong Zhang   cm = c->m;
1414b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1415e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1416b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1417b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1418b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1419b1683b59SHong Zhang     c->ncolumns[i] = n;
1420b1683b59SHong Zhang     if (n) {
1421d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1422b1683b59SHong Zhang     }
1423d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1424d2853540SHong Zhang     columns_i       += n;
1425b1683b59SHong Zhang 
1426b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1427e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1428e99be685SHong Zhang 
1429b1683b59SHong Zhang     /* loop over columns*/
1430b1683b59SHong Zhang     for (j=0; j<n; j++) {
1431b1683b59SHong Zhang       col     = is[j];
1432d2853540SHong Zhang       row_idx = cj + ci[col];
1433b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1434b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1435b1683b59SHong Zhang       for (k=0; k<m; k++) {
1436e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1437d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1438b1683b59SHong Zhang       }
1439b1683b59SHong Zhang     }
1440b1683b59SHong Zhang     /* count the number of hits */
1441b1683b59SHong Zhang     nrows = 0;
1442e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1443b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1444b1683b59SHong Zhang     }
1445b1683b59SHong Zhang     c->nrows[i]      = nrows;
1446d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1447d2853540SHong Zhang 
1448b1683b59SHong Zhang     nrows       = 0;
1449e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1450b1683b59SHong Zhang       if (rowhit[j]) {
1451d2853540SHong Zhang         rows_i[nrows]            = j;
1452e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1453b1683b59SHong Zhang         nrows++;
1454b1683b59SHong Zhang       }
1455b1683b59SHong Zhang     }
1456b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1457d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1458b1683b59SHong Zhang   }
1459e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1460b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1461b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1462d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1463d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1464d2853540SHong Zhang #endif
1465b28d80bdSHong Zhang 
1466d2853540SHong Zhang   c->colorforrow     = colorforrow;
1467d2853540SHong Zhang   c->rows            = rows;
1468d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1469d2853540SHong Zhang   c->colorforcol     = colorforcol;
1470d2853540SHong Zhang   c->columns         = columns;
1471f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1472b1683b59SHong Zhang   PetscFunctionReturn(0);
1473b1683b59SHong Zhang }
1474