xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 25c417975dded2198c01f9c84d63aaa1e09cb164)
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>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
11e9e4536cSHong Zhang /*
12e9e4536cSHong Zhang #define DEBUG_MATMATMULT
13e9e4536cSHong Zhang  */
14ec01deb9SMatthew Knepley EXTERN_C_BEGIN
156284ec50SHong Zhang #undef __FUNCT__
165c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1738baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1838baddfdSBarry Smith {
19dfbe8321SBarry Smith   PetscErrorCode ierr;
20d16893f4SBarry Smith   PetscBool      scalable=PETSC_FALSE,scalable_new=PETSC_FALSE;
215c66b693SKris Buschelman 
225c66b693SKris Buschelman   PetscFunctionBegin;
2326be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
24152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
25152983bfSHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
265e6967c8SHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable_new","Use a scalable but slower C=A*B","",scalable_new,&scalable_new,PETSC_NULL);CHKERRQ(ierr);
27d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
28d8bbc50fSBarry Smith     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2925296bd5SBarry Smith     if (scalable_new){
3025296bd5SBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_new(A,B,fill,C);CHKERRQ(ierr);
31aacf854cSBarry Smith     } else if (scalable){
32aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
3325023028SHong Zhang     } else {
3426be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3525023028SHong Zhang     }
36598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3726be0446SHong Zhang   }
38598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3901e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
40598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
415c66b693SKris Buschelman   PetscFunctionReturn(0);
425c66b693SKris Buschelman }
43ec01deb9SMatthew Knepley EXTERN_C_END
441c24bd37SHong Zhang 
45e9e4536cSHong Zhang /*
46e9e4536cSHong Zhang  MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ - Get symbolic structure of C=A*B
47e9e4536cSHong Zhang   Input Parameter:
48e9e4536cSHong Zhang .    am, Ai, Aj - number of rows and structure of A
49e9e4536cSHong Zhang .    bm, bn, Bi, Bj - number of rows, columns, and structure of B
50e9e4536cSHong Zhang .    fill - filll ratio See MatMatMult()
51e9e4536cSHong Zhang 
52e9e4536cSHong Zhang   Output Parameter:
53e9e4536cSHong Zhang .    Ci, Cj - structure of C = A*B
54e9e4536cSHong Zhang .    nspacedouble - number of extra mallocs
55e9e4536cSHong Zhang  */
5625296bd5SBarry Smith /*
5725296bd5SBarry Smith      Does not use bitarray
5825296bd5SBarry Smith  */
5925296bd5SBarry Smith #undef __FUNCT__
6025296bd5SBarry Smith #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable_new"
6125296bd5SBarry Smith PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable_new(Mat A,Mat B,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
6225296bd5SBarry Smith {
6325296bd5SBarry Smith   PetscErrorCode     ierr;
6425296bd5SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
6525296bd5SBarry Smith   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data;
6625296bd5SBarry Smith   const PetscInt     *Ai=a->i,*Aj=a->j,*Bi=b->i,*Bj=b->j,am=A->rmap->N,bm=B->rmap->N,bn=B->cmap->N;
6725296bd5SBarry Smith   const PetscInt     *aj=Aj,*bi=Bi,*bj;
6825296bd5SBarry Smith   PetscInt           *ci,*cj;
6925296bd5SBarry Smith   PetscInt           i,j,anzi,brow,bnzj,cnzi,crmax,ndouble=0;
7025296bd5SBarry Smith   PetscInt           lnk_max=bn,*lnk;
7125296bd5SBarry Smith 
7225296bd5SBarry Smith   PetscFunctionBegin;
7325296bd5SBarry Smith   /* Allocate ci array, arrays for fill computation and */
7425296bd5SBarry Smith   /* free space for accumulating nonzero column info */
7525296bd5SBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
7625296bd5SBarry Smith   ci[0] = 0;
7725296bd5SBarry Smith 
7825296bd5SBarry Smith   /* create and initialize a condensed linked list */
7925296bd5SBarry Smith   crmax = a->rmax*b->rmax;
8025296bd5SBarry Smith   lnk_max = PetscMin(crmax,bn);
8125296bd5SBarry Smith   if (!lnk_max && bn) {
8225296bd5SBarry Smith     lnk_max = bn; /* in case rmax is not defined for A or B */
8325296bd5SBarry Smith #if defined(PETSC_USE_DEBUGGING)
8425296bd5SBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: Armax or Brmax is not defined in MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable!\n");
8525296bd5SBarry Smith #endif
8625296bd5SBarry Smith   }
8725296bd5SBarry Smith #if defined(DEBUG_MATMATMULT)
8825296bd5SBarry Smith   ierr = (PETSC_COMM_SELF,"LLCondensedCreate lnk_max=%d, bn %d, crmax %d\n",lnk_max,bn,crmax);
8925296bd5SBarry Smith #endif
90aacf854cSBarry Smith   ierr = PetscLLCondensedCreate_fast(lnk_max,&lnk);CHKERRQ(ierr);
9125296bd5SBarry Smith 
9225296bd5SBarry Smith   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
9325296bd5SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(Ai[am]+Bi[bm])),&free_space);CHKERRQ(ierr);
9425296bd5SBarry Smith   current_space = free_space;
9525296bd5SBarry Smith 
9625296bd5SBarry Smith   /* Determine ci and cj */
9725296bd5SBarry Smith   for (i=0; i<am; i++) {
9825296bd5SBarry Smith     anzi = Ai[i+1] - Ai[i];
9925296bd5SBarry Smith     cnzi = 0;
10025296bd5SBarry Smith     aj   = Aj + Ai[i];
10125296bd5SBarry Smith     for (j=0; j<anzi; j++){
10225296bd5SBarry Smith       brow = aj[j];
10325296bd5SBarry Smith       bnzj = bi[brow+1] - bi[brow];
10425296bd5SBarry Smith       bj   = Bj + bi[brow];
10525296bd5SBarry Smith       /* add non-zero cols of B into the condensed sorted linked list lnk */
10625296bd5SBarry Smith       /* ierr = PetscLLCondensedView(lnk_max,lnk_max,lnk,lnk);CHKERRQ(ierr); */
107aacf854cSBarry Smith       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
10825296bd5SBarry Smith     }
109aacf854cSBarry Smith     cnzi    = lnk[1];
11025296bd5SBarry Smith     ci[i+1] = ci[i] + cnzi;
11125296bd5SBarry Smith 
11225296bd5SBarry Smith     /* If free space is not available, make more free space */
11325296bd5SBarry Smith     /* Double the amount of total space in the list */
11425296bd5SBarry Smith     if (current_space->local_remaining<cnzi) {
11525296bd5SBarry Smith       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
11625296bd5SBarry Smith       ndouble++;
11725296bd5SBarry Smith     }
11825296bd5SBarry Smith 
11925296bd5SBarry Smith     /* Copy data into free space, then initialize lnk */
120d8bbc50fSBarry Smith     /* ierr = PetscLLCondensedView_fast(lnk);CHKERRQ(ierr);   */
121aacf854cSBarry Smith     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
122aacf854cSBarry Smith     /* PetscIntView(cnzi,current_space->array,0); */
12325296bd5SBarry Smith     current_space->array           += cnzi;
12425296bd5SBarry Smith     current_space->local_used      += cnzi;
12525296bd5SBarry Smith     current_space->local_remaining -= cnzi;
12625296bd5SBarry Smith   }
12725296bd5SBarry Smith 
12825296bd5SBarry Smith   /* Column indices are in the list of free space */
12925296bd5SBarry Smith   /* Allocate and copy column indices to cj, then destroy list of free space, lnk and bt */
13025296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
13125296bd5SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
132aacf854cSBarry Smith   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
13325296bd5SBarry Smith 
13425296bd5SBarry Smith   *Ci           = ci;
13525296bd5SBarry Smith   *Cj           = cj;
13625296bd5SBarry Smith   *nspacedouble = ndouble;
13725296bd5SBarry Smith   PetscFunctionReturn(0);
13825296bd5SBarry Smith }
13925296bd5SBarry Smith 
140e9e4536cSHong Zhang #undef __FUNCT__
141b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
142b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
143b561aa9dSHong Zhang {
144b561aa9dSHong Zhang   PetscErrorCode     ierr;
145b561aa9dSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
146971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
14725c41797SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
148b561aa9dSHong Zhang   PetscReal          afill;
149fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
150fb908031SHong Zhang   PetscBT            lnkbt;
151fb908031SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
152b561aa9dSHong Zhang 
153b561aa9dSHong Zhang   PetscFunctionBegin;
154bd958071SHong Zhang   /* Get ci and cj */
155bd958071SHong Zhang   /*---------------*/
156fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
157fb908031SHong Zhang   /* free space for accumulating nonzero column info */
158fb908031SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
159fb908031SHong Zhang   ci[0] = 0;
160fb908031SHong Zhang 
161fb908031SHong Zhang   /* create and initialize a linked list */
162fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
163fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
164fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
165fb908031SHong Zhang 
166fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
167fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
168fb908031SHong Zhang   current_space = free_space;
169fb908031SHong Zhang 
170fb908031SHong Zhang   /* Determine ci and cj */
171fb908031SHong Zhang   for (i=0; i<am; i++) {
172fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
173fb908031SHong Zhang     aj   = a->j + ai[i];
174fb908031SHong Zhang     for (j=0; j<anzi; j++){
175fb908031SHong Zhang       brow = aj[j];
176fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
177fb908031SHong Zhang       bj   = b->j + bi[brow];
178fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
179fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
180fb908031SHong Zhang     }
181fb908031SHong Zhang     cnzi = lnk[0];
182fb908031SHong Zhang 
183fb908031SHong Zhang     /* If free space is not available, make more free space */
184fb908031SHong Zhang     /* Double the amount of total space in the list */
185fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
186fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
187fb908031SHong Zhang       ndouble++;
188fb908031SHong Zhang     }
189fb908031SHong Zhang 
190fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
191fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
192fb908031SHong Zhang     current_space->array           += cnzi;
193fb908031SHong Zhang     current_space->local_used      += cnzi;
194fb908031SHong Zhang     current_space->local_remaining -= cnzi;
195fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
196fb908031SHong Zhang   }
197fb908031SHong Zhang 
198fb908031SHong Zhang   /* Column indices are in the list of free space */
199fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
200fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
201fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
202fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
203fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
204b561aa9dSHong Zhang 
20526be0446SHong Zhang   /* put together the new symbolic matrix */
206aa1aec99SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr);
20758c24d83SHong Zhang 
20858c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
20958c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
21058c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
211aa1aec99SHong Zhang   c->free_a  = PETSC_FALSE;
212e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
21358c24d83SHong Zhang   c->nonew   = 0;
214002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
2150b7e3e3dSHong Zhang 
2167212da7cSHong Zhang   /* set MatInfo */
2177212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
218f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
2197212da7cSHong Zhang   c->maxnz                     = ci[am];
2207212da7cSHong Zhang   c->nz                        = ci[am];
221fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
2227212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
2237212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
2247212da7cSHong Zhang 
2257212da7cSHong Zhang #if defined(PETSC_USE_INFO)
2267212da7cSHong Zhang   if (ci[am]) {
227fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
228f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
229f2b054eeSHong Zhang   } else {
230f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
231be0fcf8dSHong Zhang   }
232f2b054eeSHong Zhang #endif
23358c24d83SHong Zhang   PetscFunctionReturn(0);
23458c24d83SHong Zhang }
235d50806bdSBarry Smith 
23626be0446SHong Zhang #undef __FUNCT__
23726be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
238dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
239d50806bdSBarry Smith {
240dfbe8321SBarry Smith   PetscErrorCode ierr;
241d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
242d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
243d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
244d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
24538baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
246c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
247519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
248aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
249aa1aec99SHong Zhang   PetscScalar    *ab_dense;
250d50806bdSBarry Smith 
251d50806bdSBarry Smith   PetscFunctionBegin;
252aa1aec99SHong Zhang   /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */
253aa1aec99SHong Zhang   if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
254aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
255aa1aec99SHong Zhang     c->a      = ca;
256aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
257aa1aec99SHong Zhang 
258aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
259aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
260aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
261aa1aec99SHong Zhang   } else {
262aa1aec99SHong Zhang     ca       = c->a;
263aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
264aa1aec99SHong Zhang   }
265aa1aec99SHong Zhang 
266c124e916SHong Zhang   /* clean old values in C */
267c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
268d50806bdSBarry Smith   /* Traverse A row-wise. */
269d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
270d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
271d50806bdSBarry Smith   for (i=0; i<am; i++) {
272d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
273d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
274519eb980SHong Zhang       brow = aj[j];
275d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
276d50806bdSBarry Smith       bjj  = bj + bi[brow];
277d50806bdSBarry Smith       baj  = ba + bi[brow];
278519eb980SHong Zhang       /* perform dense axpy */
27936ec6d2dSHong Zhang       valtmp = aa[j];
280519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
28136ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
282519eb980SHong Zhang       }
283519eb980SHong Zhang       flops += 2*bnzi;
284519eb980SHong Zhang     }
285c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
286c58ca2e3SHong Zhang 
287c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
288519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
289519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
290519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
291519eb980SHong Zhang     }
292519eb980SHong Zhang     flops += cnzi;
293519eb980SHong Zhang     cj += cnzi; ca += cnzi;
294519eb980SHong Zhang   }
295c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
296c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
297c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
298c58ca2e3SHong Zhang   PetscFunctionReturn(0);
299c58ca2e3SHong Zhang }
300c58ca2e3SHong Zhang 
301c58ca2e3SHong Zhang #undef __FUNCT__
30225023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
30325023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
304c58ca2e3SHong Zhang {
305c58ca2e3SHong Zhang   PetscErrorCode ierr;
306c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
307c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
308c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
309c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
310c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
311c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
312c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
31336ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
314c58ca2e3SHong Zhang   PetscInt       nextb;
315c58ca2e3SHong Zhang 
316c58ca2e3SHong Zhang   PetscFunctionBegin;
317e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
318971236abSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr);
319e9e4536cSHong Zhang #endif
320c58ca2e3SHong Zhang   /* clean old values in C */
321c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
322c58ca2e3SHong Zhang   /* Traverse A row-wise. */
323c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
324c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
325519eb980SHong Zhang   for (i=0;i<am;i++) {
326519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
327519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
328519eb980SHong Zhang     for (j=0;j<anzi;j++) {
329519eb980SHong Zhang       brow = aj[j];
330519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
331519eb980SHong Zhang       bjj  = bj + bi[brow];
332519eb980SHong Zhang       baj  = ba + bi[brow];
333519eb980SHong Zhang       /* perform sparse axpy */
33436ec6d2dSHong Zhang       valtmp = aa[j];
335c124e916SHong Zhang       nextb  = 0;
336c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
337c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
33836ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
339c124e916SHong Zhang         }
340d50806bdSBarry Smith       }
341d50806bdSBarry Smith       flops += 2*bnzi;
342d50806bdSBarry Smith     }
343519eb980SHong Zhang     aj += anzi; aa += anzi;
344519eb980SHong Zhang     cj += cnzi; ca += cnzi;
345d50806bdSBarry Smith   }
346c58ca2e3SHong Zhang 
347716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
348716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
349d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
350d50806bdSBarry Smith   PetscFunctionReturn(0);
351d50806bdSBarry Smith }
352bc011b1eSHong Zhang 
353e9e4536cSHong Zhang #undef __FUNCT__
35425296bd5SBarry Smith #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_new"
35525296bd5SBarry Smith PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_new(Mat A,Mat B,PetscReal fill,Mat *C)
35625296bd5SBarry Smith {
35725296bd5SBarry Smith   PetscErrorCode ierr;
35825296bd5SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
35925296bd5SBarry Smith   PetscInt       *ai=a->i,*bi=b->i,*ci,*cj;
36025296bd5SBarry Smith   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
36125296bd5SBarry Smith   MatScalar      *ca;
36225296bd5SBarry Smith   PetscReal      afill;
36325296bd5SBarry Smith 
36425296bd5SBarry Smith   PetscFunctionBegin;
36525296bd5SBarry Smith #if defined(DEBUG_MATMATMULT)
36625296bd5SBarry Smith   ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable Armax %d, Brmax %d\n",a->rmax,b->rmax);CHKERRQ(ierr);
36725296bd5SBarry Smith #endif
36825296bd5SBarry Smith   /* Get ci and cj */
36925296bd5SBarry Smith   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable_new(A,B,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
37025296bd5SBarry Smith 
37125296bd5SBarry Smith   /* Allocate space for ca */
37225296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
37325296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
37425296bd5SBarry Smith 
37525296bd5SBarry Smith   /* put together the new symbolic matrix */
37625296bd5SBarry Smith   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
37725296bd5SBarry Smith 
37825296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
37925296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
38025296bd5SBarry Smith   c = (Mat_SeqAIJ *)((*C)->data);
38125296bd5SBarry Smith   c->free_a   = PETSC_TRUE;
38225296bd5SBarry Smith   c->free_ij  = PETSC_TRUE;
38325296bd5SBarry Smith   c->nonew    = 0;
38425296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
38525296bd5SBarry Smith 
38625296bd5SBarry Smith   /* set MatInfo */
38725296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
38825296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
38925296bd5SBarry Smith   c->maxnz                     = ci[am];
39025296bd5SBarry Smith   c->nz                        = ci[am];
39125296bd5SBarry Smith   (*C)->info.mallocs           = nspacedouble;
39225296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
39325296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
39425296bd5SBarry Smith 
39525296bd5SBarry Smith #if defined(PETSC_USE_INFO)
39625296bd5SBarry Smith   if (ci[am]) {
39725296bd5SBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
39825296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
39925296bd5SBarry Smith   } else {
40025296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
40125296bd5SBarry Smith   }
40225296bd5SBarry Smith #endif
40325296bd5SBarry Smith   PetscFunctionReturn(0);
40425296bd5SBarry Smith }
40525296bd5SBarry Smith 
40625296bd5SBarry Smith 
40725296bd5SBarry Smith #undef __FUNCT__
40825023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
40925023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
410e9e4536cSHong Zhang {
411e9e4536cSHong Zhang   PetscErrorCode     ierr;
412e9e4536cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
413bf9555e6SHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
41425c41797SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
415e9e4536cSHong Zhang   MatScalar          *ca;
416e9e4536cSHong Zhang   PetscReal          afill;
417bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
418bd958071SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
419e9e4536cSHong Zhang 
420e9e4536cSHong Zhang   PetscFunctionBegin;
421bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
422bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
423bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
424bd958071SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
425bd958071SHong Zhang   ci[0] = 0;
426bd958071SHong Zhang 
427bd958071SHong Zhang   /* create and initialize a linked list */
428bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
429bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
430bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
431bd958071SHong Zhang 
432bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
433bd958071SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
434bd958071SHong Zhang   current_space = free_space;
435bd958071SHong Zhang 
436bd958071SHong Zhang   /* Determine ci and cj */
437bd958071SHong Zhang   for (i=0; i<am; i++) {
438bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
439bd958071SHong Zhang     aj   = a->j + ai[i];
440bd958071SHong Zhang     for (j=0; j<anzi; j++){
441bd958071SHong Zhang       brow = aj[j];
442bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
443bd958071SHong Zhang       bj   = b->j + bi[brow];
444bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
445bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
446bd958071SHong Zhang     }
447bd958071SHong Zhang     cnzi = lnk[0];
448bd958071SHong Zhang 
449bd958071SHong Zhang     /* If free space is not available, make more free space */
450bd958071SHong Zhang     /* Double the amount of total space in the list */
451bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
452bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
453bd958071SHong Zhang       ndouble++;
454bd958071SHong Zhang     }
455bd958071SHong Zhang 
456bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
457bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
458bd958071SHong Zhang     current_space->array           += cnzi;
459bd958071SHong Zhang     current_space->local_used      += cnzi;
460bd958071SHong Zhang     current_space->local_remaining -= cnzi;
461bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
462bd958071SHong Zhang   }
463bd958071SHong Zhang 
464bd958071SHong Zhang   /* Column indices are in the list of free space */
465bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
466bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
467bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
468bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
469bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
470e9e4536cSHong Zhang 
471e9e4536cSHong Zhang   /* Allocate space for ca */
472bd958071SHong Zhang   /*-----------------------*/
473e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
474e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
475e9e4536cSHong Zhang 
476e9e4536cSHong Zhang   /* put together the new symbolic matrix */
477e9e4536cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
478e9e4536cSHong Zhang 
479e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
480e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
481e9e4536cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
482e9e4536cSHong Zhang   c->free_a   = PETSC_TRUE;
483e9e4536cSHong Zhang   c->free_ij  = PETSC_TRUE;
484e9e4536cSHong Zhang   c->nonew    = 0;
48525023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
486e9e4536cSHong Zhang 
487e9e4536cSHong Zhang   /* set MatInfo */
488e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
489e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
490e9e4536cSHong Zhang   c->maxnz                     = ci[am];
491e9e4536cSHong Zhang   c->nz                        = ci[am];
492bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
493e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
494e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
495e9e4536cSHong Zhang 
496e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
497e9e4536cSHong Zhang   if (ci[am]) {
498bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
499e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
500e9e4536cSHong Zhang   } else {
501e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
502e9e4536cSHong Zhang   }
503e9e4536cSHong Zhang #endif
504e9e4536cSHong Zhang   PetscFunctionReturn(0);
505e9e4536cSHong Zhang }
506e9e4536cSHong Zhang 
507e9e4536cSHong Zhang 
508d2853540SHong Zhang /* This routine is not used. Should be removed! */
509bc011b1eSHong Zhang #undef __FUNCT__
5106fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
5116fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
5125df89d91SHong Zhang {
513bc011b1eSHong Zhang   PetscErrorCode ierr;
514bc011b1eSHong Zhang 
515bc011b1eSHong Zhang   PetscFunctionBegin;
516bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5176fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
518bc011b1eSHong Zhang   }
5196fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
520bc011b1eSHong Zhang   PetscFunctionReturn(0);
521bc011b1eSHong Zhang }
522bc011b1eSHong Zhang 
523bc011b1eSHong Zhang #undef __FUNCT__
524e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
525e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
5262128a86cSHong Zhang {
5272128a86cSHong Zhang   PetscErrorCode      ierr;
528e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
5292128a86cSHong Zhang 
5302128a86cSHong Zhang   PetscFunctionBegin;
531b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
5322128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
5332128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
5342128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
5352128a86cSHong Zhang   PetscFunctionReturn(0);
5362128a86cSHong Zhang }
5372128a86cSHong Zhang 
5382128a86cSHong Zhang #undef __FUNCT__
5392128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
5402128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
5412128a86cSHong Zhang {
5422128a86cSHong Zhang   PetscErrorCode      ierr;
5432128a86cSHong Zhang   PetscContainer      container;
544e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
5452128a86cSHong Zhang 
5462128a86cSHong Zhang   PetscFunctionBegin;
547e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
5482128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
5492128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
5502128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
5512128a86cSHong Zhang   if (A->ops->destroy) {
5522128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
5532128a86cSHong Zhang   }
554e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
5552128a86cSHong Zhang   PetscFunctionReturn(0);
5562128a86cSHong Zhang }
5572128a86cSHong Zhang 
5582128a86cSHong Zhang #undef __FUNCT__
5596fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
5606fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
561bc011b1eSHong Zhang {
5625df89d91SHong Zhang   PetscErrorCode      ierr;
56381d82fe4SHong Zhang   Mat                 Bt;
56481d82fe4SHong Zhang   PetscInt            *bti,*btj;
565e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
5662128a86cSHong Zhang   PetscContainer      container;
567d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
568d2853540SHong Zhang 
56981d82fe4SHong Zhang   PetscFunctionBegin;
570d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
57181d82fe4SHong Zhang    /* create symbolic Bt */
57281d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
57381d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
57481d82fe4SHong Zhang 
57581d82fe4SHong Zhang   /* get symbolic C=A*Bt */
57681d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
57781d82fe4SHong Zhang 
5782128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
579e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
5802128a86cSHong Zhang 
5812128a86cSHong Zhang   /* attach the supporting struct to C */
5822128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
5832128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
584e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
585e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
5862128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5872128a86cSHong Zhang 
5882128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
5892128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
5902128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
5912128a86cSHong Zhang 
592d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
593d2853540SHong Zhang   etime2 += tf - t0;
594d2853540SHong Zhang 
595f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
5962128a86cSHong Zhang   if (multtrans->usecoloring){
597b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
598b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
5992128a86cSHong Zhang     ISColoring           iscoloring;
6002128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
601d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
602e8354b3eSHong Zhang 
603e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
604502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
605d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
606d2853540SHong Zhang     etime0 += tf - t0;
607d2853540SHong Zhang 
608d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
609b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
6102128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
6112128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
612e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
613d2853540SHong Zhang     etime01 += tf - t0;
6142128a86cSHong Zhang 
615e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
6162128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
6172128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
6182128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
6192128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
6202128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
6212128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
6222128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
6232128a86cSHong Zhang 
6242128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
6252128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
6262128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
6272128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
6282128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
6292128a86cSHong Zhang     multtrans->ABt_den = C_dense;
630e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
631e8354b3eSHong Zhang     etime1 += tf - t0;
632f94ccd6cSHong Zhang 
633f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
6341ce71dffSSatish Balay     {
635f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
636f94ccd6cSHong 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));
637f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
6381ce71dffSSatish Balay     }
639f94ccd6cSHong Zhang #endif
6402128a86cSHong Zhang   }
641*e99be685SHong Zhang   /* clean up */
642*e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
643*e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
6442128a86cSHong Zhang 
645f94ccd6cSHong Zhang 
646f94ccd6cSHong Zhang 
64781d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
64881d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
6491fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
6501fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6511fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
6521fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
6531fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
6541fa1209cSHong Zhang   MatScalar          *ca;
6551fa1209cSHong Zhang   PetscBT            lnkbt;
6561fa1209cSHong Zhang   PetscReal          afill;
6575df89d91SHong Zhang 
6581fa1209cSHong Zhang   /* Allocate row pointer array ci  */
6591fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6601fa1209cSHong Zhang   ci[0] = 0;
6611fa1209cSHong Zhang 
6621fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
6631fa1209cSHong Zhang   nlnk = bm+1;
6641fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
6651fa1209cSHong Zhang 
6661fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
6671fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
6681fa1209cSHong Zhang   current_space = free_space;
6691fa1209cSHong Zhang 
6701fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
6711fa1209cSHong Zhang   for (i=0; i<am; i++) {
6721fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
6731fa1209cSHong Zhang     cnzi = 0;
6741fa1209cSHong Zhang     acol = aj + ai[i];
6751fa1209cSHong Zhang     for (j=0; j<bm; j++){
6761fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
6771fa1209cSHong Zhang       bcol= bj + bi[j];
67881d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
6791fa1209cSHong Zhang       ka = 0; kb = 0;
6801fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
68181d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
68281d82fe4SHong Zhang         if (ka == anzi) break;
68381d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
68481d82fe4SHong Zhang         if (kb == bnzj) break;
6851fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
6861fa1209cSHong Zhang           index[0] = j;
6871fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
6881fa1209cSHong Zhang           cnzi++;
6891fa1209cSHong Zhang           break;
6901fa1209cSHong Zhang         }
6911fa1209cSHong Zhang       }
6921fa1209cSHong Zhang     }
6931fa1209cSHong Zhang 
6941fa1209cSHong Zhang     /* If free space is not available, make more free space */
6951fa1209cSHong Zhang     /* Double the amount of total space in the list */
6961fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
6971fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
6981fa1209cSHong Zhang       nspacedouble++;
6991fa1209cSHong Zhang     }
7001fa1209cSHong Zhang 
7011fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
7021fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
7031fa1209cSHong Zhang     current_space->array           += cnzi;
7041fa1209cSHong Zhang     current_space->local_used      += cnzi;
7051fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
7061fa1209cSHong Zhang 
7071fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
7081fa1209cSHong Zhang   }
7091fa1209cSHong Zhang 
7101fa1209cSHong Zhang 
7111fa1209cSHong Zhang   /* Column indices are in the list of free space.
7121fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
7131fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
7141fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
7151fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7161fa1209cSHong Zhang 
7171fa1209cSHong Zhang   /* Allocate space for ca */
7181fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
7191fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
7201fa1209cSHong Zhang 
7211fa1209cSHong Zhang   /* put together the new symbolic matrix */
7221fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
7231fa1209cSHong Zhang 
7241fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
7251fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
7261fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
7271fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
7281fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
7291fa1209cSHong Zhang   c->nonew    = 0;
7301fa1209cSHong Zhang 
7311fa1209cSHong Zhang   /* set MatInfo */
7321fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
7331fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
7341fa1209cSHong Zhang   c->maxnz                     = ci[am];
7351fa1209cSHong Zhang   c->nz                        = ci[am];
7361fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
7371fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
7381fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
7391fa1209cSHong Zhang 
7401fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
7411fa1209cSHong Zhang   if (ci[am]) {
7421fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
7436fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
7441fa1209cSHong Zhang   } else {
7451fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7461fa1209cSHong Zhang   }
7471fa1209cSHong Zhang #endif
74881d82fe4SHong Zhang #endif
7495df89d91SHong Zhang   PetscFunctionReturn(0);
7505df89d91SHong Zhang }
7515df89d91SHong Zhang 
7526973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
7535df89d91SHong Zhang #undef __FUNCT__
7546fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
7556fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
7565df89d91SHong Zhang {
7575df89d91SHong Zhang   PetscErrorCode ierr;
7585df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
759e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
76081d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
7615df89d91SHong Zhang   PetscLogDouble flops=0.0;
762aa1aec99SHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval;
763e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
7642128a86cSHong Zhang   PetscContainer      container;
7656973769bSHong Zhang #if defined(USE_ARRAY)
7666973769bSHong Zhang   MatScalar      *spdot;
7676973769bSHong Zhang #endif
7685df89d91SHong Zhang 
7695df89d91SHong Zhang   PetscFunctionBegin;
770e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
7712128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
7722128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
7732128a86cSHong Zhang   if (multtrans->usecoloring){
774b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
775c8db22f6SHong Zhang     Mat                   Bt_dense;
776c8db22f6SHong Zhang     PetscInt              m,n;
777b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
778a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
779a0b95ffbSSatish Balay 
7802128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
781c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
782c8db22f6SHong Zhang 
783b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
7842128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
785b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
7862128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
787b2d2b04fSHong Zhang     etime0 += tf - t0;
788c8db22f6SHong Zhang 
789c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
7902128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
7912128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
7922128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7932128a86cSHong Zhang     etime2 += tf - t0;
794c8db22f6SHong Zhang 
7952128a86cSHong Zhang     /* Recover C from C_dense */
7962128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
797b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
7982128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
7992128a86cSHong Zhang     etime1 += tf - t0;
800f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
801f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
802f94ccd6cSHong Zhang #endif
803c8db22f6SHong Zhang     PetscFunctionReturn(0);
804c8db22f6SHong Zhang   }
805c8db22f6SHong Zhang 
8066973769bSHong Zhang #if defined(USE_ARRAY)
8076973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
808e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
809e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
8106973769bSHong Zhang #endif
8116973769bSHong Zhang 
81281d82fe4SHong Zhang   /* clear old values in C */
813aa1aec99SHong Zhang   if (!c->a){
814aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
815aa1aec99SHong Zhang     c->a      = ca;
816aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
817aa1aec99SHong Zhang   } else {
818aa1aec99SHong Zhang     ca =  c->a;
819aa1aec99SHong Zhang   }
82081d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
8215df89d91SHong Zhang 
8221fa1209cSHong Zhang   for (i=0; i<cm; i++) {
82381d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
82481d82fe4SHong Zhang     acol = aj + ai[i];
8256973769bSHong Zhang     aval = aa + ai[i];
8261fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
8271fa1209cSHong Zhang     ccol = cj + ci[i];
8286973769bSHong Zhang     cval = ca + ci[i];
8291fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
83081d82fe4SHong Zhang       brow = ccol[j];
83181d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
83281d82fe4SHong Zhang       bcol = bj + bi[brow];
8336973769bSHong Zhang       bval = ba + bi[brow];
8346973769bSHong Zhang 
83581d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
8366973769bSHong Zhang #if defined(USE_ARRAY)
8376973769bSHong Zhang       /* put ba to spdot array */
8386973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
8396973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
8406973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
8416973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
8426973769bSHong Zhang       }
8436973769bSHong Zhang       /* zero spdot array */
8446973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
8456973769bSHong Zhang #else
84681d82fe4SHong Zhang       nexta = 0; nextb = 0;
84781d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
84881d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
84981d82fe4SHong Zhang         if (nexta == anzi) break;
85081d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
85181d82fe4SHong Zhang         if (nextb == bnzj) break;
85281d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
8536973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
85481d82fe4SHong Zhang           nexta++; nextb++;
85581d82fe4SHong Zhang           flops += 2;
8561fa1209cSHong Zhang         }
8571fa1209cSHong Zhang       }
8586973769bSHong Zhang #endif
85981d82fe4SHong Zhang     }
86081d82fe4SHong Zhang   }
86181d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
86281d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
86381d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
8646973769bSHong Zhang #if defined(USE_ARRAY)
8656973769bSHong Zhang   ierr = PetscFree(spdot);
8666973769bSHong Zhang #endif
8675df89d91SHong Zhang   PetscFunctionReturn(0);
8685df89d91SHong Zhang }
8695df89d91SHong Zhang 
8705df89d91SHong Zhang #undef __FUNCT__
87175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
87275648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
8735df89d91SHong Zhang   PetscErrorCode ierr;
8745df89d91SHong Zhang 
8755df89d91SHong Zhang   PetscFunctionBegin;
8765df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
87775648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8785df89d91SHong Zhang   }
87975648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8805df89d91SHong Zhang   PetscFunctionReturn(0);
8815df89d91SHong Zhang }
8825df89d91SHong Zhang 
8835df89d91SHong Zhang #undef __FUNCT__
88475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
88575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
8865df89d91SHong Zhang {
887bc011b1eSHong Zhang   PetscErrorCode ierr;
888bc011b1eSHong Zhang   Mat            At;
88938baddfdSBarry Smith   PetscInt       *ati,*atj;
890bc011b1eSHong Zhang 
891bc011b1eSHong Zhang   PetscFunctionBegin;
892bc011b1eSHong Zhang   /* create symbolic At */
893bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
894d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
895bc011b1eSHong Zhang 
896bc011b1eSHong Zhang   /* get symbolic C=At*B */
897bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
898bc011b1eSHong Zhang 
899bc011b1eSHong Zhang   /* clean up */
9006bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
901bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
902bc011b1eSHong Zhang   PetscFunctionReturn(0);
903bc011b1eSHong Zhang }
904bc011b1eSHong Zhang 
905bc011b1eSHong Zhang #undef __FUNCT__
90675648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
90775648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
908bc011b1eSHong Zhang {
909bc011b1eSHong Zhang   PetscErrorCode ierr;
9100fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
911d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
912d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
913d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
914aa1aec99SHong Zhang   MatScalar      *aa=a->a,*ba,*ca,*caj;
915bc011b1eSHong Zhang 
916bc011b1eSHong Zhang   PetscFunctionBegin;
917aa1aec99SHong Zhang   if (!c->a){
918aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
919aa1aec99SHong Zhang     c->a      = ca;
920aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
921aa1aec99SHong Zhang   } else {
922aa1aec99SHong Zhang     ca = c->a;
923aa1aec99SHong Zhang   }
924bc011b1eSHong Zhang   /* clear old values in C */
925bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
926bc011b1eSHong Zhang 
927bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
928bc011b1eSHong Zhang   for (i=0;i<am;i++) {
929bc011b1eSHong Zhang     bj   = b->j + bi[i];
930bc011b1eSHong Zhang     ba   = b->a + bi[i];
931bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
932bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
933bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
934bc011b1eSHong Zhang       nextb = 0;
9350fbc74f4SHong Zhang       crow  = *aj++;
936bc011b1eSHong Zhang       cjj   = cj + ci[crow];
937bc011b1eSHong Zhang       caj   = ca + ci[crow];
938bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
939bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
9400fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
9410fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
942bc011b1eSHong Zhang           nextb++;
943bc011b1eSHong Zhang         }
944bc011b1eSHong Zhang       }
945bc011b1eSHong Zhang       flops += 2*bnzi;
9460fbc74f4SHong Zhang       aa++;
947bc011b1eSHong Zhang     }
948bc011b1eSHong Zhang   }
949bc011b1eSHong Zhang 
950bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
951bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
952bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
953bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
954bc011b1eSHong Zhang   PetscFunctionReturn(0);
955bc011b1eSHong Zhang }
9569123193aSHong Zhang 
957ec01deb9SMatthew Knepley EXTERN_C_BEGIN
9589123193aSHong Zhang #undef __FUNCT__
9599123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
9609123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9619123193aSHong Zhang {
9629123193aSHong Zhang   PetscErrorCode ierr;
9639123193aSHong Zhang 
9649123193aSHong Zhang   PetscFunctionBegin;
9659123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
9669123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
9679123193aSHong Zhang   }
9689123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
9699123193aSHong Zhang   PetscFunctionReturn(0);
9709123193aSHong Zhang }
971ec01deb9SMatthew Knepley EXTERN_C_END
972edd81eecSMatthew Knepley 
9739123193aSHong Zhang #undef __FUNCT__
9749123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
9759123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
9769123193aSHong Zhang {
9779123193aSHong Zhang   PetscErrorCode ierr;
9789123193aSHong Zhang 
9799123193aSHong Zhang   PetscFunctionBegin;
9805a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
9818cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
9829123193aSHong Zhang   PetscFunctionReturn(0);
9839123193aSHong Zhang }
9849123193aSHong Zhang 
9859123193aSHong Zhang #undef __FUNCT__
9869123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
9879123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
9889123193aSHong Zhang {
989f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
9909123193aSHong Zhang   PetscErrorCode ierr;
991dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
992dd6ea824SBarry Smith   MatScalar      *aa;
993d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
994f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
9959123193aSHong Zhang 
9969123193aSHong Zhang   PetscFunctionBegin;
997f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
998e32f2f54SBarry 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);
999e32f2f54SBarry 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);
1000e32f2f54SBarry 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);
1001f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
1002f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
1003f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1004f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
1005f32d5d43SBarry Smith     colam = col*am;
1006f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1007f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1008f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
1009f32d5d43SBarry Smith       aj  = a->j + a->i[i];
1010f32d5d43SBarry Smith       aa  = a->a + a->i[i];
1011f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1012f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
1013f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
1014f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
1015f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
10169123193aSHong Zhang       }
1017f32d5d43SBarry Smith       c[colam + i]       = r1;
1018f32d5d43SBarry Smith       c[colam + am + i]  = r2;
1019f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
1020f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
1021f32d5d43SBarry Smith     }
1022f32d5d43SBarry Smith     b1 += bm4;
1023f32d5d43SBarry Smith     b2 += bm4;
1024f32d5d43SBarry Smith     b3 += bm4;
1025f32d5d43SBarry Smith     b4 += bm4;
1026f32d5d43SBarry Smith   }
1027f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
1028f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1029f32d5d43SBarry Smith       r1 = 0.0;
1030f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
1031f32d5d43SBarry Smith       aj  = a->j + a->i[i];
1032f32d5d43SBarry Smith       aa  = a->a + a->i[i];
1033f32d5d43SBarry Smith 
1034f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1035f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
1036f32d5d43SBarry Smith       }
1037f32d5d43SBarry Smith       c[col*am + i]     = r1;
1038f32d5d43SBarry Smith     }
1039f32d5d43SBarry Smith     b1 += bm;
1040f32d5d43SBarry Smith   }
1041dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
1042f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1043f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
1044f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1045f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1046f32d5d43SBarry Smith   PetscFunctionReturn(0);
1047f32d5d43SBarry Smith }
1048f32d5d43SBarry Smith 
1049f32d5d43SBarry Smith /*
10504324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1051f32d5d43SBarry Smith */
1052f32d5d43SBarry Smith #undef __FUNCT__
1053f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1054f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1055f32d5d43SBarry Smith {
1056f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1057f32d5d43SBarry Smith   PetscErrorCode ierr;
1058dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1059dd6ea824SBarry Smith   MatScalar      *aa;
1060d0f46423SBarry 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;
10614324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1062f32d5d43SBarry Smith 
1063f32d5d43SBarry Smith   PetscFunctionBegin;
1064f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1065f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
1066f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
1067f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
10684324174eSBarry Smith 
10694324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
10704324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
10714324174eSBarry Smith       colam = col*am;
10724324174eSBarry Smith       arm   = a->compressedrow.nrows;
10734324174eSBarry Smith       ii    = a->compressedrow.i;
10744324174eSBarry Smith       ridx  = a->compressedrow.rindex;
10754324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
10764324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
10774324174eSBarry Smith 	n   = ii[i+1] - ii[i];
10784324174eSBarry Smith 	aj  = a->j + ii[i];
10794324174eSBarry Smith 	aa  = a->a + ii[i];
10804324174eSBarry Smith 	for (j=0; j<n; j++) {
10814324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
10824324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
10834324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
10844324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
10854324174eSBarry Smith 	}
10864324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
10874324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
10884324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
10894324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
10904324174eSBarry Smith       }
10914324174eSBarry Smith       b1 += bm4;
10924324174eSBarry Smith       b2 += bm4;
10934324174eSBarry Smith       b3 += bm4;
10944324174eSBarry Smith       b4 += bm4;
10954324174eSBarry Smith     }
10964324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
10974324174eSBarry Smith       colam = col*am;
10984324174eSBarry Smith       arm   = a->compressedrow.nrows;
10994324174eSBarry Smith       ii    = a->compressedrow.i;
11004324174eSBarry Smith       ridx  = a->compressedrow.rindex;
11014324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
11024324174eSBarry Smith 	r1 = 0.0;
11034324174eSBarry Smith 	n   = ii[i+1] - ii[i];
11044324174eSBarry Smith 	aj  = a->j + ii[i];
11054324174eSBarry Smith 	aa  = a->a + ii[i];
11064324174eSBarry Smith 
11074324174eSBarry Smith 	for (j=0; j<n; j++) {
11084324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
11094324174eSBarry Smith 	}
11104324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
11114324174eSBarry Smith       }
11124324174eSBarry Smith       b1 += bm;
11134324174eSBarry Smith     }
11144324174eSBarry Smith   } else {
1115f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
1116f32d5d43SBarry Smith       colam = col*am;
1117f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1118f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
1119f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1120f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1121f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1122f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1123f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
1124f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
1125f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
1126f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
1127f32d5d43SBarry Smith 	}
1128f32d5d43SBarry Smith 	c[colam + i]       += r1;
1129f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
1130f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
1131f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
1132f32d5d43SBarry Smith       }
1133f32d5d43SBarry Smith       b1 += bm4;
1134f32d5d43SBarry Smith       b2 += bm4;
1135f32d5d43SBarry Smith       b3 += bm4;
1136f32d5d43SBarry Smith       b4 += bm4;
1137f32d5d43SBarry Smith     }
1138f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
1139f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1140f32d5d43SBarry Smith 	r1 = 0.0;
1141f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1142f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1143f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1144f32d5d43SBarry Smith 
1145f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1146f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
1147f32d5d43SBarry Smith 	}
1148f32d5d43SBarry Smith 	c[col*am + i]     += r1;
1149f32d5d43SBarry Smith       }
1150f32d5d43SBarry Smith       b1 += bm;
1151f32d5d43SBarry Smith     }
11524324174eSBarry Smith   }
1153dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
1154f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1155f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
11569123193aSHong Zhang   PetscFunctionReturn(0);
11579123193aSHong Zhang }
1158b1683b59SHong Zhang 
1159b1683b59SHong Zhang #undef __FUNCT__
1160b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1161b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1162c8db22f6SHong Zhang {
1163c8db22f6SHong Zhang   PetscErrorCode ierr;
11642f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
11652f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
11662f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
11672f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
11682f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
11692f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1170c8db22f6SHong Zhang 
1171c8db22f6SHong Zhang   PetscFunctionBegin;
11722f5f1f90SHong Zhang   btval_den=btdense->v;
11732f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
11742f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
11752f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
11762f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1177d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
11782f5f1f90SHong Zhang       btcol = bj + bi[col];
11792f5f1f90SHong Zhang       btval = ba + bi[col];
11802f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1181d2853540SHong Zhang       for (j=0; j<anz; j++){
11822f5f1f90SHong Zhang         brow            = btcol[j];
11832f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1184c8db22f6SHong Zhang       }
1185c8db22f6SHong Zhang     }
11862f5f1f90SHong Zhang     btval_den += m;
1187c8db22f6SHong Zhang   }
1188c8db22f6SHong Zhang   PetscFunctionReturn(0);
1189c8db22f6SHong Zhang }
1190c8db22f6SHong Zhang 
1191c8db22f6SHong Zhang #undef __FUNCT__
1192b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1193b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
11948972f759SHong Zhang {
11958972f759SHong Zhang   PetscErrorCode ierr;
1196b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
11972f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1198b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
11992f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
12008972f759SHong Zhang 
12018972f759SHong Zhang   PetscFunctionBegin;
12028972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
1203b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
1204b2d2b04fSHong Zhang   cp_den = ca_den;
12052f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
12062f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
12072f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
12082f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
12092f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
12102f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1211b2d2b04fSHong Zhang     }
1212b2d2b04fSHong Zhang     cp_den += m;
1213b2d2b04fSHong Zhang   }
1214b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
12158972f759SHong Zhang   PetscFunctionReturn(0);
12168972f759SHong Zhang }
12178972f759SHong Zhang 
1218*e99be685SHong Zhang /*
1219*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1220*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1221*e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1222*e99be685SHong Zhang  */
1223*e99be685SHong Zhang #undef __FUNCT__
1224*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
1225*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1226*e99be685SHong Zhang {
1227*e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1228*e99be685SHong Zhang   PetscErrorCode ierr;
1229*e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1230*e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1231*e99be685SHong Zhang   PetscInt       *cspidx;
1232*e99be685SHong Zhang 
1233*e99be685SHong Zhang   PetscFunctionBegin;
1234*e99be685SHong Zhang   *nn = n;
1235*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1236*e99be685SHong Zhang   if (symmetric) {
1237*e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
1238*e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
1239*e99be685SHong Zhang   } else {
1240*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1241*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1242*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1243*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1244*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1245*e99be685SHong Zhang     jj = a->j;
1246*e99be685SHong Zhang     for (i=0; i<nz; i++) {
1247*e99be685SHong Zhang       collengths[jj[i]]++;
1248*e99be685SHong Zhang     }
1249*e99be685SHong Zhang     cia[0] = oshift;
1250*e99be685SHong Zhang     for (i=0; i<n; i++) {
1251*e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1252*e99be685SHong Zhang     }
1253*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1254*e99be685SHong Zhang     jj   = a->j;
1255*e99be685SHong Zhang     for (row=0; row<m; row++) {
1256*e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1257*e99be685SHong Zhang       for (i=0; i<mr; i++) {
1258*e99be685SHong Zhang         col = *jj++;
1259*e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1260*e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1261*e99be685SHong Zhang       }
1262*e99be685SHong Zhang     }
1263*e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1264*e99be685SHong Zhang     *ia = cia; *ja = cja;
1265*e99be685SHong Zhang     *spidx = cspidx;
1266*e99be685SHong Zhang   }
1267*e99be685SHong Zhang   PetscFunctionReturn(0);
1268*e99be685SHong Zhang }
1269*e99be685SHong Zhang 
1270*e99be685SHong Zhang #undef __FUNCT__
1271*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1272*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1273*e99be685SHong Zhang {
1274*e99be685SHong Zhang   PetscErrorCode ierr;
1275*e99be685SHong Zhang 
1276*e99be685SHong Zhang   PetscFunctionBegin;
1277*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1278*e99be685SHong Zhang 
1279*e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1280*e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1281*e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1282*e99be685SHong Zhang   PetscFunctionReturn(0);
1283*e99be685SHong Zhang }
1284*e99be685SHong Zhang 
12858972f759SHong Zhang #undef __FUNCT__
1286b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1287b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1288b1683b59SHong Zhang {
1289b1683b59SHong Zhang   PetscErrorCode ierr;
1290d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1291b1683b59SHong Zhang   const PetscInt *is;
1292b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1293b1683b59SHong Zhang   IS             *isa;
1294b1683b59SHong Zhang   PetscBool      done;
1295b1683b59SHong Zhang   PetscBool      flg1,flg2;
1296b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1297*e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1298d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1299b1683b59SHong Zhang 
1300b1683b59SHong Zhang   PetscFunctionBegin;
1301b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1302*e99be685SHong Zhang 
1303b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1304b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1305b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1306b1683b59SHong Zhang   if (flg1 || flg2) {
1307b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1308b1683b59SHong Zhang   }
1309b1683b59SHong Zhang 
1310b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1311b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1312b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1313b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1314b1683b59SHong Zhang   c->rstart = 0;
1315b1683b59SHong Zhang 
1316b1683b59SHong Zhang   c->ncolors = nis;
1317b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1318b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1319d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1320d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1321d2853540SHong Zhang   colorforrow[0]    = 0;
1322d2853540SHong Zhang   rows_i            = rows;
1323d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1324b1683b59SHong Zhang 
1325d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1326d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1327d2853540SHong Zhang   colorforcol[0] = 0;
1328d2853540SHong Zhang   columns_i      = columns;
1329d2853540SHong Zhang 
1330*e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1331b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1332b1683b59SHong Zhang 
1333b28d80bdSHong Zhang   cm = c->m;
1334b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1335*e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1336b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1337b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1338b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1339b1683b59SHong Zhang     c->ncolumns[i] = n;
1340b1683b59SHong Zhang     if (n) {
1341d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1342b1683b59SHong Zhang     }
1343d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1344d2853540SHong Zhang     columns_i       += n;
1345b1683b59SHong Zhang 
1346b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1347e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1348*e99be685SHong Zhang 
1349b1683b59SHong Zhang     /* loop over columns*/
1350b1683b59SHong Zhang     for (j=0; j<n; j++) {
1351b1683b59SHong Zhang       col     = is[j];
1352d2853540SHong Zhang       row_idx = cj + ci[col];
1353b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1354b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1355b1683b59SHong Zhang       for (k=0; k<m; k++) {
1356*e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1357d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1358b1683b59SHong Zhang       }
1359b1683b59SHong Zhang     }
1360b1683b59SHong Zhang     /* count the number of hits */
1361b1683b59SHong Zhang     nrows = 0;
1362e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1363b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1364b1683b59SHong Zhang     }
1365b1683b59SHong Zhang     c->nrows[i]      = nrows;
1366d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1367d2853540SHong Zhang 
1368b1683b59SHong Zhang     nrows       = 0;
1369e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1370b1683b59SHong Zhang       if (rowhit[j]) {
1371d2853540SHong Zhang         rows_i[nrows]            = j;
1372*e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1373b1683b59SHong Zhang         nrows++;
1374b1683b59SHong Zhang       }
1375b1683b59SHong Zhang     }
1376b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1377d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1378b1683b59SHong Zhang   }
1379*e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1380b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1381b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1382d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1383d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1384d2853540SHong Zhang #endif
1385b28d80bdSHong Zhang 
1386d2853540SHong Zhang   c->colorforrow     = colorforrow;
1387d2853540SHong Zhang   c->rows            = rows;
1388d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1389d2853540SHong Zhang   c->colorforcol     = colorforcol;
1390d2853540SHong Zhang   c->columns         = columns;
1391f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1392b1683b59SHong Zhang   PetscFunctionReturn(0);
1393b1683b59SHong Zhang }
1394