xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 07475bc16356fc37e8c66fcce1957fb7d8feef24)
1be1d678aSKris Buschelman 
2d50806bdSBarry Smith /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4d50806bdSBarry Smith           C = A * B
5d50806bdSBarry Smith */
6d50806bdSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
90ced3a2bSJed Brown #include <../src/mat/utils/petscheap.h>
10c6db04a5SJed Brown #include <petscbt.h>
11*07475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h>
127bab7c10SHong Zhang 
136284ec50SHong Zhang #undef __FUNCT__
145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1638baddfdSBarry Smith {
17dfbe8321SBarry Smith   PetscErrorCode ierr;
188a07c6f1SJed Brown   PetscBool      scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE,btheap = PETSC_FALSE;
195c66b693SKris Buschelman 
205c66b693SKris Buschelman   PetscFunctionBegin;
2126be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
22152983bfSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
23152983bfSHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
243c50cad2SHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,PETSC_NULL);CHKERRQ(ierr);
250ced3a2bSJed Brown     ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,PETSC_NULL);CHKERRQ(ierr);
268a07c6f1SJed Brown     ierr = PetscOptionsBool("-matmatmult_btheap","Use btheap implementation of symbolic factorization C=A*B","",btheap,&btheap,PETSC_NULL);CHKERRQ(ierr);
27d8bbc50fSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
283c50cad2SHong Zhang     if (scalable_fast){
293c50cad2SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
30aacf854cSBarry Smith     } else if (scalable){
31aacf854cSBarry Smith       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
320ced3a2bSJed Brown     } else if (heap) {
330ced3a2bSJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
348a07c6f1SJed Brown     } else if (btheap) {
358a07c6f1SJed Brown       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
3625023028SHong Zhang     } else {
3726be0446SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3825023028SHong Zhang     }
3926be0446SHong Zhang   }
405c913ed7SHong Zhang 
4101e47db4SHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
425c66b693SKris Buschelman   PetscFunctionReturn(0);
435c66b693SKris Buschelman }
441c24bd37SHong Zhang 
45e9e4536cSHong Zhang #undef __FUNCT__
46b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
47b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
48b561aa9dSHong Zhang {
49b561aa9dSHong Zhang   PetscErrorCode     ierr;
50b561aa9dSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
51971236abSHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
52c1ba5575SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
53b561aa9dSHong Zhang   PetscReal          afill;
54fb908031SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
55fb908031SHong Zhang   PetscBT            lnkbt;
56fb908031SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
57b561aa9dSHong Zhang 
58b561aa9dSHong Zhang   PetscFunctionBegin;
59bd958071SHong Zhang   /* Get ci and cj */
60bd958071SHong Zhang   /*---------------*/
61fb908031SHong Zhang   /* Allocate ci array, arrays for fill computation and */
62fb908031SHong Zhang   /* free space for accumulating nonzero column info */
63fb908031SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
64fb908031SHong Zhang   ci[0] = 0;
65fb908031SHong Zhang 
66fb908031SHong Zhang   /* create and initialize a linked list */
67fb908031SHong Zhang   nlnk_max = a->rmax*b->rmax;
68fb908031SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
69fb908031SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr);
70fb908031SHong Zhang 
71fb908031SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
72fb908031SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
73fb908031SHong Zhang   current_space = free_space;
74fb908031SHong Zhang 
75fb908031SHong Zhang   /* Determine ci and cj */
76fb908031SHong Zhang   for (i=0; i<am; i++) {
77fb908031SHong Zhang     anzi = ai[i+1] - ai[i];
78fb908031SHong Zhang     aj   = a->j + ai[i];
79fb908031SHong Zhang     for (j=0; j<anzi; j++){
80fb908031SHong Zhang       brow = aj[j];
81fb908031SHong Zhang       bnzj = bi[brow+1] - bi[brow];
82fb908031SHong Zhang       bj   = b->j + bi[brow];
83fb908031SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
84fb908031SHong Zhang       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
85fb908031SHong Zhang     }
86fb908031SHong Zhang     cnzi = lnk[0];
87fb908031SHong Zhang 
88fb908031SHong Zhang     /* If free space is not available, make more free space */
89fb908031SHong Zhang     /* Double the amount of total space in the list */
90fb908031SHong Zhang     if (current_space->local_remaining<cnzi) {
91fb908031SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
92fb908031SHong Zhang       ndouble++;
93fb908031SHong Zhang     }
94fb908031SHong Zhang 
95fb908031SHong Zhang     /* Copy data into free space, then initialize lnk */
96fb908031SHong Zhang     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
97fb908031SHong Zhang     current_space->array           += cnzi;
98fb908031SHong Zhang     current_space->local_used      += cnzi;
99fb908031SHong Zhang     current_space->local_remaining -= cnzi;
100fb908031SHong Zhang     ci[i+1] = ci[i] + cnzi;
101fb908031SHong Zhang   }
102fb908031SHong Zhang 
103fb908031SHong Zhang   /* Column indices are in the list of free space */
104fb908031SHong Zhang   /* Allocate space for cj, initialize cj, and */
105fb908031SHong Zhang   /* destroy list of free space and other temporary array(s) */
106fb908031SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
107fb908031SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
108fb908031SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
109b561aa9dSHong Zhang 
11026be0446SHong Zhang   /* put together the new symbolic matrix */
111aa1aec99SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr);
112a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
113a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
11458c24d83SHong Zhang 
11558c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
11658c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
11758c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
118aa1aec99SHong Zhang   c->free_a  = PETSC_FALSE;
119e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
12058c24d83SHong Zhang   c->nonew   = 0;
121002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1220b7e3e3dSHong Zhang 
1237212da7cSHong Zhang   /* set MatInfo */
1247212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
125f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1267212da7cSHong Zhang   c->maxnz                     = ci[am];
1277212da7cSHong Zhang   c->nz                        = ci[am];
128fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1297212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1307212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1317212da7cSHong Zhang 
1327212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1337212da7cSHong Zhang   if (ci[am]) {
134fb908031SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
135f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
136f2b054eeSHong Zhang   } else {
137f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
138be0fcf8dSHong Zhang   }
139f2b054eeSHong Zhang #endif
14058c24d83SHong Zhang   PetscFunctionReturn(0);
14158c24d83SHong Zhang }
142d50806bdSBarry Smith 
14326be0446SHong Zhang #undef __FUNCT__
14426be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
145dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
146d50806bdSBarry Smith {
147dfbe8321SBarry Smith   PetscErrorCode ierr;
148d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
149d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
150d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
151d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
15238baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
153c58ca2e3SHong Zhang   PetscInt       am=A->rmap->n,cm=C->rmap->n;
154519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
155aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
156aa1aec99SHong Zhang   PetscScalar    *ab_dense;
157d50806bdSBarry Smith 
158d50806bdSBarry Smith   PetscFunctionBegin;
159aa1aec99SHong Zhang   /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */
160aa1aec99SHong Zhang   if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
161aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
162aa1aec99SHong Zhang     c->a      = ca;
163aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
164aa1aec99SHong Zhang 
165aa1aec99SHong Zhang     ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr);
166aa1aec99SHong Zhang     ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
167aa1aec99SHong Zhang     c->matmult_abdense = ab_dense;
168aa1aec99SHong Zhang   } else {
169aa1aec99SHong Zhang     ca       = c->a;
170aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
171aa1aec99SHong Zhang   }
172aa1aec99SHong Zhang 
173c124e916SHong Zhang   /* clean old values in C */
174c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
175d50806bdSBarry Smith   /* Traverse A row-wise. */
176d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
177d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
178d50806bdSBarry Smith   for (i=0; i<am; i++) {
179d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
180d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
181519eb980SHong Zhang       brow = aj[j];
182d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
183d50806bdSBarry Smith       bjj  = bj + bi[brow];
184d50806bdSBarry Smith       baj  = ba + bi[brow];
185519eb980SHong Zhang       /* perform dense axpy */
18636ec6d2dSHong Zhang       valtmp = aa[j];
187519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
18836ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
189519eb980SHong Zhang       }
190519eb980SHong Zhang       flops += 2*bnzi;
191519eb980SHong Zhang     }
192c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
193c58ca2e3SHong Zhang 
194c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
195519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
196519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
197519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
198519eb980SHong Zhang     }
199519eb980SHong Zhang     flops += cnzi;
200519eb980SHong Zhang     cj += cnzi; ca += cnzi;
201519eb980SHong Zhang   }
202c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
203c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
204c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
205c58ca2e3SHong Zhang   PetscFunctionReturn(0);
206c58ca2e3SHong Zhang }
207c58ca2e3SHong Zhang 
208c58ca2e3SHong Zhang #undef __FUNCT__
20925023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable"
21025023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
211c58ca2e3SHong Zhang {
212c58ca2e3SHong Zhang   PetscErrorCode ierr;
213c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
214c58ca2e3SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
215c58ca2e3SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
216c58ca2e3SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
217c58ca2e3SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
218c58ca2e3SHong Zhang   PetscInt       am=A->rmap->N,cm=C->rmap->N;
219c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
22036ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
221c58ca2e3SHong Zhang   PetscInt       nextb;
222c58ca2e3SHong Zhang 
223c58ca2e3SHong Zhang   PetscFunctionBegin;
224c58ca2e3SHong Zhang   /* clean old values in C */
225c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
226c58ca2e3SHong Zhang   /* Traverse A row-wise. */
227c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
228c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
229519eb980SHong Zhang   for (i=0;i<am;i++) {
230519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
231519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
232519eb980SHong Zhang     for (j=0;j<anzi;j++) {
233519eb980SHong Zhang       brow = aj[j];
234519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
235519eb980SHong Zhang       bjj  = bj + bi[brow];
236519eb980SHong Zhang       baj  = ba + bi[brow];
237519eb980SHong Zhang       /* perform sparse axpy */
23836ec6d2dSHong Zhang       valtmp = aa[j];
239c124e916SHong Zhang       nextb  = 0;
240c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
241c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
24236ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
243c124e916SHong Zhang         }
244d50806bdSBarry Smith       }
245d50806bdSBarry Smith       flops += 2*bnzi;
246d50806bdSBarry Smith     }
247519eb980SHong Zhang     aj += anzi; aa += anzi;
248519eb980SHong Zhang     cj += cnzi; ca += cnzi;
249d50806bdSBarry Smith   }
250c58ca2e3SHong Zhang 
251716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
254d50806bdSBarry Smith   PetscFunctionReturn(0);
255d50806bdSBarry Smith }
256bc011b1eSHong Zhang 
257e9e4536cSHong Zhang #undef __FUNCT__
2583c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast"
2593c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
26025296bd5SBarry Smith {
26125296bd5SBarry Smith   PetscErrorCode     ierr;
26225296bd5SBarry Smith   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
26325296bd5SBarry Smith   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
2643c50cad2SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
26525296bd5SBarry Smith   MatScalar          *ca;
26625296bd5SBarry Smith   PetscReal          afill;
2673c50cad2SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
2683c50cad2SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
26925296bd5SBarry Smith 
27025296bd5SBarry Smith   PetscFunctionBegin;
2713c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
2723c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
2733c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
2743c50cad2SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2753c50cad2SHong Zhang   ci[0] = 0;
2763c50cad2SHong Zhang 
2773c50cad2SHong Zhang   /* create and initialize a linked list */
2783c50cad2SHong Zhang   nlnk_max = a->rmax*b->rmax;
2793c50cad2SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
2803c50cad2SHong Zhang   ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr);
2813c50cad2SHong Zhang 
2823c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
2833c50cad2SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
2843c50cad2SHong Zhang   current_space = free_space;
2853c50cad2SHong Zhang 
2863c50cad2SHong Zhang   /* Determine ci and cj */
2873c50cad2SHong Zhang   for (i=0; i<am; i++) {
2883c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
2893c50cad2SHong Zhang     aj   = a->j + ai[i];
2903c50cad2SHong Zhang     for (j=0; j<anzi; j++){
2913c50cad2SHong Zhang       brow = aj[j];
2923c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
2933c50cad2SHong Zhang       bj   = b->j + bi[brow];
2943c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
2953c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
2963c50cad2SHong Zhang     }
2973c50cad2SHong Zhang     cnzi = lnk[1];
2983c50cad2SHong Zhang 
2993c50cad2SHong Zhang     /* If free space is not available, make more free space */
3003c50cad2SHong Zhang     /* Double the amount of total space in the list */
3013c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
3023c50cad2SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3033c50cad2SHong Zhang       ndouble++;
3043c50cad2SHong Zhang     }
3053c50cad2SHong Zhang 
3063c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3073c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3083c50cad2SHong Zhang     current_space->array           += cnzi;
3093c50cad2SHong Zhang     current_space->local_used      += cnzi;
3103c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3113c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3123c50cad2SHong Zhang   }
3133c50cad2SHong Zhang 
3143c50cad2SHong Zhang   /* Column indices are in the list of free space */
3153c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3163c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
3173c50cad2SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3183c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3193c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
32025296bd5SBarry Smith 
32125296bd5SBarry Smith   /* Allocate space for ca */
32225296bd5SBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
32325296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
32425296bd5SBarry Smith 
32525296bd5SBarry Smith   /* put together the new symbolic matrix */
32625296bd5SBarry Smith   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
327a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
328a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
32925296bd5SBarry Smith 
33025296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
33125296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
33225296bd5SBarry Smith   c = (Mat_SeqAIJ *)((*C)->data);
33325296bd5SBarry Smith   c->free_a   = PETSC_TRUE;
33425296bd5SBarry Smith   c->free_ij  = PETSC_TRUE;
33525296bd5SBarry Smith   c->nonew    = 0;
33625296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
33725296bd5SBarry Smith 
33825296bd5SBarry Smith   /* set MatInfo */
33925296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
34025296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
34125296bd5SBarry Smith   c->maxnz                     = ci[am];
34225296bd5SBarry Smith   c->nz                        = ci[am];
3433c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
34425296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
34525296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
34625296bd5SBarry Smith 
34725296bd5SBarry Smith #if defined(PETSC_USE_INFO)
34825296bd5SBarry Smith   if (ci[am]) {
3493c50cad2SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
35025296bd5SBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
35125296bd5SBarry Smith   } else {
35225296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
35325296bd5SBarry Smith   }
35425296bd5SBarry Smith #endif
35525296bd5SBarry Smith   PetscFunctionReturn(0);
35625296bd5SBarry Smith }
35725296bd5SBarry Smith 
35825296bd5SBarry Smith 
35925296bd5SBarry Smith #undef __FUNCT__
36025023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable"
36125023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
362e9e4536cSHong Zhang {
363e9e4536cSHong Zhang   PetscErrorCode     ierr;
364e9e4536cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
365bf9555e6SHong Zhang   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
36625c41797SHong Zhang   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
367e9e4536cSHong Zhang   MatScalar          *ca;
368e9e4536cSHong Zhang   PetscReal          afill;
369bd958071SHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
370bd958071SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
371e9e4536cSHong Zhang 
372e9e4536cSHong Zhang   PetscFunctionBegin;
373bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
374bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
375bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
376bd958071SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
377bd958071SHong Zhang   ci[0] = 0;
378bd958071SHong Zhang 
379bd958071SHong Zhang   /* create and initialize a linked list */
380bd958071SHong Zhang   nlnk_max = a->rmax*b->rmax;
381bd958071SHong Zhang   if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
382bd958071SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
383bd958071SHong Zhang 
384bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
385bd958071SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
386bd958071SHong Zhang   current_space = free_space;
387bd958071SHong Zhang 
388bd958071SHong Zhang   /* Determine ci and cj */
389bd958071SHong Zhang   for (i=0; i<am; i++) {
390bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
391bd958071SHong Zhang     aj   = a->j + ai[i];
392bd958071SHong Zhang     for (j=0; j<anzi; j++){
393bd958071SHong Zhang       brow = aj[j];
394bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
395bd958071SHong Zhang       bj   = b->j + bi[brow];
396bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
397bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
398bd958071SHong Zhang     }
399bd958071SHong Zhang     cnzi = lnk[0];
400bd958071SHong Zhang 
401bd958071SHong Zhang     /* If free space is not available, make more free space */
402bd958071SHong Zhang     /* Double the amount of total space in the list */
403bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
404bd958071SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
405bd958071SHong Zhang       ndouble++;
406bd958071SHong Zhang     }
407bd958071SHong Zhang 
408bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
409bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
410bd958071SHong Zhang     current_space->array           += cnzi;
411bd958071SHong Zhang     current_space->local_used      += cnzi;
412bd958071SHong Zhang     current_space->local_remaining -= cnzi;
413bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
414bd958071SHong Zhang   }
415bd958071SHong Zhang 
416bd958071SHong Zhang   /* Column indices are in the list of free space */
417bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
418bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
419bd958071SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
420bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
421bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
422e9e4536cSHong Zhang 
423e9e4536cSHong Zhang   /* Allocate space for ca */
424bd958071SHong Zhang   /*-----------------------*/
425e9e4536cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
426e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
427e9e4536cSHong Zhang 
428e9e4536cSHong Zhang   /* put together the new symbolic matrix */
429e9e4536cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
430a2f3521dSMark F. Adams   (*C)->rmap->bs = A->rmap->bs;
431a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
432e9e4536cSHong Zhang 
433e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
434e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
435e9e4536cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
436e9e4536cSHong Zhang   c->free_a   = PETSC_TRUE;
437e9e4536cSHong Zhang   c->free_ij  = PETSC_TRUE;
438e9e4536cSHong Zhang   c->nonew    = 0;
43925023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
440e9e4536cSHong Zhang 
441e9e4536cSHong Zhang   /* set MatInfo */
442e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
443e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
444e9e4536cSHong Zhang   c->maxnz                     = ci[am];
445e9e4536cSHong Zhang   c->nz                        = ci[am];
446bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
447e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
448e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
449e9e4536cSHong Zhang 
450e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
451e9e4536cSHong Zhang   if (ci[am]) {
452bd958071SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
453e9e4536cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
454e9e4536cSHong Zhang   } else {
455e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
456e9e4536cSHong Zhang   }
457e9e4536cSHong Zhang #endif
458e9e4536cSHong Zhang   PetscFunctionReturn(0);
459e9e4536cSHong Zhang }
460e9e4536cSHong Zhang 
4610ced3a2bSJed Brown #undef __FUNCT__
4620ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap"
4630ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
4640ced3a2bSJed Brown {
4650ced3a2bSJed Brown   PetscErrorCode     ierr;
4660ced3a2bSJed Brown   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
4670ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
4680ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
4690ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
4700ced3a2bSJed Brown   PetscReal          afill;
4710ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
4720ced3a2bSJed Brown   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
4730ced3a2bSJed Brown   PetscHeap          h;
4740ced3a2bSJed Brown 
4750ced3a2bSJed Brown   PetscFunctionBegin;
4760ced3a2bSJed Brown   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
4770ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
4780ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
4790ced3a2bSJed Brown   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4800ced3a2bSJed Brown   ci[0] = 0;
4810ced3a2bSJed Brown 
4820ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
4830ced3a2bSJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
4840ced3a2bSJed Brown   current_space = free_space;
4850ced3a2bSJed Brown 
4860ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
4870ced3a2bSJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
4880ced3a2bSJed Brown 
4890ced3a2bSJed Brown   /* Determine ci and cj */
4900ced3a2bSJed Brown   for (i=0; i<am; i++) {
4910ced3a2bSJed 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 */
4920ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
4930ced3a2bSJed Brown     ci[i+1] = ci[i];
4940ced3a2bSJed Brown     /* Populate the min heap */
4950ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
4960ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
4970ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
4980ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
4990ced3a2bSJed Brown       }
5000ced3a2bSJed Brown     }
5010ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5020ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5030ced3a2bSJed Brown     while (j >= 0) {
5040ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
5050ced3a2bSJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
5060ced3a2bSJed Brown         ndouble++;
5070ced3a2bSJed Brown       }
5080ced3a2bSJed Brown       *(current_space->array++) = col;
5090ced3a2bSJed Brown       current_space->local_used++;
5100ced3a2bSJed Brown       current_space->local_remaining--;
5110ced3a2bSJed Brown       ci[i+1]++;
5120ced3a2bSJed Brown 
5130ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5140ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5150ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5160ced3a2bSJed Brown         PetscInt j2,col2;
5170ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5180ced3a2bSJed Brown         if (col2 != col) break;
5190ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5200ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5210ced3a2bSJed Brown       }
5220ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5230ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5240ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5250ced3a2bSJed Brown     }
5260ced3a2bSJed Brown   }
5270ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5280ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5290ced3a2bSJed Brown 
5300ced3a2bSJed Brown   /* Column indices are in the list of free space */
5310ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5320ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
5330ced3a2bSJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5340ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5350ced3a2bSJed Brown 
5360ced3a2bSJed Brown   /* put together the new symbolic matrix */
53789d95c1aSJed Brown   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr);
5380ced3a2bSJed Brown   (*C)->rmap->bs = A->rmap->bs;
5390ced3a2bSJed Brown   (*C)->cmap->bs = B->cmap->bs;
5400ced3a2bSJed Brown 
5410ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5420ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5430ced3a2bSJed Brown   c = (Mat_SeqAIJ *)((*C)->data);
5440ced3a2bSJed Brown   c->free_a   = PETSC_TRUE;
5450ced3a2bSJed Brown   c->free_ij  = PETSC_TRUE;
5460ced3a2bSJed Brown   c->nonew    = 0;
54789d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5480ced3a2bSJed Brown 
5490ced3a2bSJed Brown   /* set MatInfo */
5500ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
5510ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
5520ced3a2bSJed Brown   c->maxnz                     = ci[am];
5530ced3a2bSJed Brown   c->nz                        = ci[am];
5540ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
5550ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
5560ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
5570ced3a2bSJed Brown 
5580ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
5590ced3a2bSJed Brown   if (ci[am]) {
5600ced3a2bSJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
5610ced3a2bSJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
5620ced3a2bSJed Brown   } else {
5630ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
5640ced3a2bSJed Brown   }
5650ced3a2bSJed Brown #endif
5660ced3a2bSJed Brown   PetscFunctionReturn(0);
5670ced3a2bSJed Brown }
568e9e4536cSHong Zhang 
5698a07c6f1SJed Brown #undef __FUNCT__
5708a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap"
5718a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
5728a07c6f1SJed Brown {
5738a07c6f1SJed Brown   PetscErrorCode     ierr;
5748a07c6f1SJed Brown   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5758a07c6f1SJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
5768a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
5778a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
5788a07c6f1SJed Brown   PetscReal          afill;
5798a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
5808a07c6f1SJed Brown   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
5818a07c6f1SJed Brown   PetscHeap          h;
5828a07c6f1SJed Brown   PetscBT            bt;
5838a07c6f1SJed Brown 
5848a07c6f1SJed Brown   PetscFunctionBegin;
5858a07c6f1SJed Brown   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
5868a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
5878a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
5888a07c6f1SJed Brown   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5898a07c6f1SJed Brown   ci[0] = 0;
5908a07c6f1SJed Brown 
5918a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
5928a07c6f1SJed Brown   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5938a07c6f1SJed Brown   current_space = free_space;
5948a07c6f1SJed Brown 
5958a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
5968a07c6f1SJed Brown   ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr);
5978a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
5988a07c6f1SJed Brown 
5998a07c6f1SJed Brown   /* Determine ci and cj */
6008a07c6f1SJed Brown   for (i=0; i<am; i++) {
6018a07c6f1SJed 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 */
6028a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6038a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6048a07c6f1SJed Brown     ci[i+1] = ci[i];
6058a07c6f1SJed Brown     /* Populate the min heap */
6068a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6078a07c6f1SJed Brown       PetscInt brow = acol[j];
6088a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6098a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6108a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6118a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6128a07c6f1SJed Brown           bb[j]++;
6138a07c6f1SJed Brown           break;
6148a07c6f1SJed Brown         }
6158a07c6f1SJed Brown       }
6168a07c6f1SJed Brown     }
6178a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6188a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6198a07c6f1SJed Brown     while (j >= 0) {
6208a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6218a07c6f1SJed Brown         fptr = PETSC_NULL;                      /* need PetscBTMemzero */
6228a07c6f1SJed Brown         ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),&current_space);CHKERRQ(ierr);
6238a07c6f1SJed Brown         ndouble++;
6248a07c6f1SJed Brown       }
6258a07c6f1SJed Brown       *(current_space->array++) = col;
6268a07c6f1SJed Brown       current_space->local_used++;
6278a07c6f1SJed Brown       current_space->local_remaining--;
6288a07c6f1SJed Brown       ci[i+1]++;
6298a07c6f1SJed Brown 
6308a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6318a07c6f1SJed Brown       for ( ; bb[j] < bi[acol[j]+1]; bb[j]++) {
6328a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6338a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6348a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6358a07c6f1SJed Brown           bb[j]++;
6368a07c6f1SJed Brown           break;
6378a07c6f1SJed Brown         }
6388a07c6f1SJed Brown       }
6398a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6408a07c6f1SJed Brown     }
6418a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6428a07c6f1SJed Brown       for ( ; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6438a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6448a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6458a07c6f1SJed Brown     }
6468a07c6f1SJed Brown   }
6478a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6488a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6498a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6508a07c6f1SJed Brown 
6518a07c6f1SJed Brown   /* Column indices are in the list of free space */
6528a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
6538a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
6548a07c6f1SJed Brown   ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6558a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
6568a07c6f1SJed Brown 
6578a07c6f1SJed Brown   /* put together the new symbolic matrix */
65889d95c1aSJed Brown   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr);
6598a07c6f1SJed Brown   (*C)->rmap->bs = A->rmap->bs;
6608a07c6f1SJed Brown   (*C)->cmap->bs = B->cmap->bs;
6618a07c6f1SJed Brown 
6628a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6638a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
6648a07c6f1SJed Brown   c = (Mat_SeqAIJ *)((*C)->data);
6658a07c6f1SJed Brown   c->free_a   = PETSC_TRUE;
6668a07c6f1SJed Brown   c->free_ij  = PETSC_TRUE;
6678a07c6f1SJed Brown   c->nonew    = 0;
66889d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
6698a07c6f1SJed Brown 
6708a07c6f1SJed Brown   /* set MatInfo */
6718a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6728a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
6738a07c6f1SJed Brown   c->maxnz                     = ci[am];
6748a07c6f1SJed Brown   c->nz                        = ci[am];
6758a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
6768a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
6778a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
6788a07c6f1SJed Brown 
6798a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
6808a07c6f1SJed Brown   if (ci[am]) {
6818a07c6f1SJed Brown     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr);
6828a07c6f1SJed Brown     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
6838a07c6f1SJed Brown   } else {
6848a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6858a07c6f1SJed Brown   }
6868a07c6f1SJed Brown #endif
6878a07c6f1SJed Brown   PetscFunctionReturn(0);
6888a07c6f1SJed Brown }
6898a07c6f1SJed Brown 
690d2853540SHong Zhang /* This routine is not used. Should be removed! */
691bc011b1eSHong Zhang #undef __FUNCT__
6926fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
6936fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
6945df89d91SHong Zhang {
695bc011b1eSHong Zhang   PetscErrorCode ierr;
696bc011b1eSHong Zhang 
697bc011b1eSHong Zhang   PetscFunctionBegin;
698bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6996fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
700bc011b1eSHong Zhang   }
7016fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
702bc011b1eSHong Zhang   PetscFunctionReturn(0);
703bc011b1eSHong Zhang }
704bc011b1eSHong Zhang 
705bc011b1eSHong Zhang #undef __FUNCT__
706e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
707e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
7082128a86cSHong Zhang {
7092128a86cSHong Zhang   PetscErrorCode      ierr;
710e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
7112128a86cSHong Zhang 
7122128a86cSHong Zhang   PetscFunctionBegin;
713b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
7142128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
7152128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
7162128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
7172128a86cSHong Zhang   PetscFunctionReturn(0);
7182128a86cSHong Zhang }
7192128a86cSHong Zhang 
7202128a86cSHong Zhang #undef __FUNCT__
7212128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
7222128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
7232128a86cSHong Zhang {
7242128a86cSHong Zhang   PetscErrorCode      ierr;
7252128a86cSHong Zhang   PetscContainer      container;
726e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
7272128a86cSHong Zhang 
7282128a86cSHong Zhang   PetscFunctionBegin;
729e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
7302128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
7312128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
7322128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
7332128a86cSHong Zhang   if (A->ops->destroy) {
7342128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
7352128a86cSHong Zhang   }
736e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
7372128a86cSHong Zhang   PetscFunctionReturn(0);
7382128a86cSHong Zhang }
7392128a86cSHong Zhang 
7402128a86cSHong Zhang #undef __FUNCT__
7416fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
7426fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
743bc011b1eSHong Zhang {
7445df89d91SHong Zhang   PetscErrorCode      ierr;
74581d82fe4SHong Zhang   Mat                 Bt;
74681d82fe4SHong Zhang   PetscInt            *bti,*btj;
747e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
7482128a86cSHong Zhang   PetscContainer      container;
749d2853540SHong Zhang 
75081d82fe4SHong Zhang   PetscFunctionBegin;
75181d82fe4SHong Zhang    /* create symbolic Bt */
75281d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
75381d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
754c3e5699bSSatish Balay   Bt->rmap->bs = A->cmap->bs;
755c3e5699bSSatish Balay   Bt->cmap->bs = B->cmap->bs;
75681d82fe4SHong Zhang 
75781d82fe4SHong Zhang   /* get symbolic C=A*Bt */
75881d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
75981d82fe4SHong Zhang 
7602128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
761e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
7622128a86cSHong Zhang 
7632128a86cSHong Zhang   /* attach the supporting struct to C */
7642128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
7652128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
766e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
767e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
7682128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
7692128a86cSHong Zhang 
7702128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
7712128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
7722128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
7732128a86cSHong Zhang 
774f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
7752128a86cSHong Zhang   if (multtrans->usecoloring){
776b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
777b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
7782128a86cSHong Zhang     ISColoring           iscoloring;
7792128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
780e8354b3eSHong Zhang 
781502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
782b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
7832128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
7842128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
7852128a86cSHong Zhang 
7862128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
7872128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
7882128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
7892128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
7902128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
7912128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
7922128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
7932128a86cSHong Zhang 
7942128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
7952128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
7962128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
7972128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
7982128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
7992128a86cSHong Zhang     multtrans->ABt_den = C_dense;
800f94ccd6cSHong Zhang 
801f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
8021ce71dffSSatish Balay     {
803f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
80422d28d08SBarry Smith     ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr);
8051ce71dffSSatish Balay     }
806f94ccd6cSHong Zhang #endif
8072128a86cSHong Zhang   }
808e99be685SHong Zhang   /* clean up */
809e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
810e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8112128a86cSHong Zhang 
812f94ccd6cSHong Zhang 
813f94ccd6cSHong Zhang 
81481d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
81581d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
8161fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
8171fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
8181fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
8191fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
8201fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
8211fa1209cSHong Zhang   MatScalar          *ca;
8221fa1209cSHong Zhang   PetscBT            lnkbt;
8231fa1209cSHong Zhang   PetscReal          afill;
8245df89d91SHong Zhang 
8251fa1209cSHong Zhang   /* Allocate row pointer array ci  */
8261fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
8271fa1209cSHong Zhang   ci[0] = 0;
8281fa1209cSHong Zhang 
8291fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
8301fa1209cSHong Zhang   nlnk = bm+1;
8311fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
8321fa1209cSHong Zhang 
8331fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
8341fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
8351fa1209cSHong Zhang   current_space = free_space;
8361fa1209cSHong Zhang 
8371fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
8381fa1209cSHong Zhang   for (i=0; i<am; i++) {
8391fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
8401fa1209cSHong Zhang     cnzi = 0;
8411fa1209cSHong Zhang     acol = aj + ai[i];
8421fa1209cSHong Zhang     for (j=0; j<bm; j++){
8431fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
8441fa1209cSHong Zhang       bcol= bj + bi[j];
84581d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
8461fa1209cSHong Zhang       ka = 0; kb = 0;
8471fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
84881d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
84981d82fe4SHong Zhang         if (ka == anzi) break;
85081d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
85181d82fe4SHong Zhang         if (kb == bnzj) break;
8521fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
8531fa1209cSHong Zhang           index[0] = j;
8541fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
8551fa1209cSHong Zhang           cnzi++;
8561fa1209cSHong Zhang           break;
8571fa1209cSHong Zhang         }
8581fa1209cSHong Zhang       }
8591fa1209cSHong Zhang     }
8601fa1209cSHong Zhang 
8611fa1209cSHong Zhang     /* If free space is not available, make more free space */
8621fa1209cSHong Zhang     /* Double the amount of total space in the list */
8631fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
8641fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
8651fa1209cSHong Zhang       nspacedouble++;
8661fa1209cSHong Zhang     }
8671fa1209cSHong Zhang 
8681fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
8691fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
8701fa1209cSHong Zhang     current_space->array           += cnzi;
8711fa1209cSHong Zhang     current_space->local_used      += cnzi;
8721fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
8731fa1209cSHong Zhang 
8741fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
8751fa1209cSHong Zhang   }
8761fa1209cSHong Zhang 
8771fa1209cSHong Zhang 
8781fa1209cSHong Zhang   /* Column indices are in the list of free space.
8791fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
8801fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
8811fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
8821fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
8831fa1209cSHong Zhang 
8841fa1209cSHong Zhang   /* Allocate space for ca */
8851fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
8861fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
8871fa1209cSHong Zhang 
8881fa1209cSHong Zhang   /* put together the new symbolic matrix */
8891fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
890a2f3521dSMark F. Adams   (*C)->rmap->bs = A->cmap->bs;
891a2f3521dSMark F. Adams   (*C)->cmap->bs = B->cmap->bs;
8921fa1209cSHong Zhang 
8931fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
8941fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
8951fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
8961fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
8971fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
8981fa1209cSHong Zhang   c->nonew    = 0;
8991fa1209cSHong Zhang 
9001fa1209cSHong Zhang   /* set MatInfo */
9011fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
9021fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
9031fa1209cSHong Zhang   c->maxnz                     = ci[am];
9041fa1209cSHong Zhang   c->nz                        = ci[am];
9051fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
9061fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
9071fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
9081fa1209cSHong Zhang 
9091fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
9101fa1209cSHong Zhang   if (ci[am]) {
9111fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
9126fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
9131fa1209cSHong Zhang   } else {
9141fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
9151fa1209cSHong Zhang   }
9161fa1209cSHong Zhang #endif
91781d82fe4SHong Zhang #endif
9185df89d91SHong Zhang   PetscFunctionReturn(0);
9195df89d91SHong Zhang }
9205df89d91SHong Zhang 
9216973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
9225df89d91SHong Zhang #undef __FUNCT__
9236fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
9246fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9255df89d91SHong Zhang {
9265df89d91SHong Zhang   PetscErrorCode ierr;
9275df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
928e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
92981d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9305df89d91SHong Zhang   PetscLogDouble flops=0.0;
931aa1aec99SHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval;
932e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
9332128a86cSHong Zhang   PetscContainer      container;
9346973769bSHong Zhang #if defined(USE_ARRAY)
9356973769bSHong Zhang   MatScalar      *spdot;
9366973769bSHong Zhang #endif
9375df89d91SHong Zhang 
9385df89d91SHong Zhang   PetscFunctionBegin;
93958ed3ceaSHong Zhang   /* clear old values in C */
94058ed3ceaSHong Zhang   if (!c->a){
94158ed3ceaSHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
94258ed3ceaSHong Zhang     c->a      = ca;
94358ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
94458ed3ceaSHong Zhang   } else {
94558ed3ceaSHong Zhang     ca =  c->a;
94658ed3ceaSHong Zhang   }
94758ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
94858ed3ceaSHong Zhang 
949e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
9502128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
9512128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
9522128a86cSHong Zhang   if (multtrans->usecoloring){
953b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
954c8db22f6SHong Zhang     Mat                   Bt_dense;
955c8db22f6SHong Zhang     PetscInt              m,n;
956a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
957a0b95ffbSSatish Balay 
9582128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
959c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
960c8db22f6SHong Zhang 
961b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
962b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
963c8db22f6SHong Zhang 
964c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9652128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
966c8db22f6SHong Zhang 
9672128a86cSHong Zhang     /* Recover C from C_dense */
968b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
969c8db22f6SHong Zhang     PetscFunctionReturn(0);
970c8db22f6SHong Zhang   }
971c8db22f6SHong Zhang 
9726973769bSHong Zhang #if defined(USE_ARRAY)
9736973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
974e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
975e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
9766973769bSHong Zhang #endif
9776973769bSHong Zhang 
9781fa1209cSHong Zhang   for (i=0; i<cm; i++) {
97981d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
98081d82fe4SHong Zhang     acol = aj + ai[i];
9816973769bSHong Zhang     aval = aa + ai[i];
9821fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9831fa1209cSHong Zhang     ccol = cj + ci[i];
9846973769bSHong Zhang     cval = ca + ci[i];
9851fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
98681d82fe4SHong Zhang       brow = ccol[j];
98781d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
98881d82fe4SHong Zhang       bcol = bj + bi[brow];
9896973769bSHong Zhang       bval = ba + bi[brow];
9906973769bSHong Zhang 
99181d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
9926973769bSHong Zhang #if defined(USE_ARRAY)
9936973769bSHong Zhang       /* put ba to spdot array */
9946973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
9956973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
9966973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
9976973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
9986973769bSHong Zhang       }
9996973769bSHong Zhang       /* zero spdot array */
10006973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
10016973769bSHong Zhang #else
100281d82fe4SHong Zhang       nexta = 0; nextb = 0;
100381d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
100481d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
100581d82fe4SHong Zhang         if (nexta == anzi) break;
100681d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
100781d82fe4SHong Zhang         if (nextb == bnzj) break;
100881d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
10096973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
101081d82fe4SHong Zhang           nexta++; nextb++;
101181d82fe4SHong Zhang           flops += 2;
10121fa1209cSHong Zhang         }
10131fa1209cSHong Zhang       }
10146973769bSHong Zhang #endif
101581d82fe4SHong Zhang     }
101681d82fe4SHong Zhang   }
101781d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101881d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101981d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10206973769bSHong Zhang #if defined(USE_ARRAY)
10216973769bSHong Zhang   ierr = PetscFree(spdot);
10226973769bSHong Zhang #endif
10235df89d91SHong Zhang   PetscFunctionReturn(0);
10245df89d91SHong Zhang }
10255df89d91SHong Zhang 
10265df89d91SHong Zhang #undef __FUNCT__
102775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
102875648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
10295df89d91SHong Zhang   PetscErrorCode ierr;
10305df89d91SHong Zhang 
10315df89d91SHong Zhang   PetscFunctionBegin;
10325df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
103375648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
10345df89d91SHong Zhang   }
103575648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
10365df89d91SHong Zhang   PetscFunctionReturn(0);
10375df89d91SHong Zhang }
10385df89d91SHong Zhang 
10395df89d91SHong Zhang #undef __FUNCT__
104075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
104175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10425df89d91SHong Zhang {
1043bc011b1eSHong Zhang   PetscErrorCode ierr;
1044bc011b1eSHong Zhang   Mat            At;
104538baddfdSBarry Smith   PetscInt       *ati,*atj;
1046bc011b1eSHong Zhang 
1047bc011b1eSHong Zhang   PetscFunctionBegin;
1048bc011b1eSHong Zhang   /* create symbolic At */
1049bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1050d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
1051a2f3521dSMark F. Adams   At->rmap->bs = A->cmap->bs;
1052a2f3521dSMark F. Adams   At->cmap->bs = B->cmap->bs;
1053bc011b1eSHong Zhang 
1054bc011b1eSHong Zhang   /* get symbolic C=At*B */
1055bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1056bc011b1eSHong Zhang 
1057bc011b1eSHong Zhang   /* clean up */
10586bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1059bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1060bc011b1eSHong Zhang   PetscFunctionReturn(0);
1061bc011b1eSHong Zhang }
1062bc011b1eSHong Zhang 
1063bc011b1eSHong Zhang #undef __FUNCT__
106475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
106575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1066bc011b1eSHong Zhang {
1067bc011b1eSHong Zhang   PetscErrorCode ierr;
10680fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1069d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1070d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1071d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1072aa1aec99SHong Zhang   MatScalar      *aa=a->a,*ba,*ca,*caj;
1073bc011b1eSHong Zhang 
1074bc011b1eSHong Zhang   PetscFunctionBegin;
1075aa1aec99SHong Zhang   if (!c->a){
1076aa1aec99SHong Zhang     ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
1077aa1aec99SHong Zhang     c->a      = ca;
1078aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1079aa1aec99SHong Zhang   } else {
1080aa1aec99SHong Zhang     ca = c->a;
1081aa1aec99SHong Zhang   }
1082bc011b1eSHong Zhang   /* clear old values in C */
1083bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1084bc011b1eSHong Zhang 
1085bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1086bc011b1eSHong Zhang   for (i=0;i<am;i++) {
1087bc011b1eSHong Zhang     bj   = b->j + bi[i];
1088bc011b1eSHong Zhang     ba   = b->a + bi[i];
1089bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1090bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1091bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1092bc011b1eSHong Zhang       nextb = 0;
10930fbc74f4SHong Zhang       crow  = *aj++;
1094bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1095bc011b1eSHong Zhang       caj   = ca + ci[crow];
1096bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1097bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
10980fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
10990fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1100bc011b1eSHong Zhang           nextb++;
1101bc011b1eSHong Zhang         }
1102bc011b1eSHong Zhang       }
1103bc011b1eSHong Zhang       flops += 2*bnzi;
11040fbc74f4SHong Zhang       aa++;
1105bc011b1eSHong Zhang     }
1106bc011b1eSHong Zhang   }
1107bc011b1eSHong Zhang 
1108bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1109bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1110bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1111bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1112bc011b1eSHong Zhang   PetscFunctionReturn(0);
1113bc011b1eSHong Zhang }
11149123193aSHong Zhang 
1115ec01deb9SMatthew Knepley EXTERN_C_BEGIN
11169123193aSHong Zhang #undef __FUNCT__
11179123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
11189123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11199123193aSHong Zhang {
11209123193aSHong Zhang   PetscErrorCode ierr;
11219123193aSHong Zhang 
11229123193aSHong Zhang   PetscFunctionBegin;
11239123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
11249123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11259123193aSHong Zhang   }
11269123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11279123193aSHong Zhang   PetscFunctionReturn(0);
11289123193aSHong Zhang }
1129ec01deb9SMatthew Knepley EXTERN_C_END
1130edd81eecSMatthew Knepley 
11319123193aSHong Zhang #undef __FUNCT__
11329123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
11339123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11349123193aSHong Zhang {
11359123193aSHong Zhang   PetscErrorCode ierr;
11369123193aSHong Zhang 
11379123193aSHong Zhang   PetscFunctionBegin;
11385a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
1139d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11409123193aSHong Zhang   PetscFunctionReturn(0);
11419123193aSHong Zhang }
11429123193aSHong Zhang 
11439123193aSHong Zhang #undef __FUNCT__
11449123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
11459123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11469123193aSHong Zhang {
1147f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11489123193aSHong Zhang   PetscErrorCode ierr;
1149dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1150dd6ea824SBarry Smith   MatScalar      *aa;
1151d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
1152f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
11539123193aSHong Zhang 
11549123193aSHong Zhang   PetscFunctionBegin;
1155f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1156e32f2f54SBarry 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);
1157e32f2f54SBarry 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);
1158e32f2f54SBarry 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);
11598c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
11608c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1161f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1162f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
1163f32d5d43SBarry Smith     colam = col*am;
1164f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1165f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1166f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
1167f32d5d43SBarry Smith       aj  = a->j + a->i[i];
1168f32d5d43SBarry Smith       aa  = a->a + a->i[i];
1169f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1170f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
1171f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
1172f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
1173f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
11749123193aSHong Zhang       }
1175f32d5d43SBarry Smith       c[colam + i]       = r1;
1176f32d5d43SBarry Smith       c[colam + am + i]  = r2;
1177f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
1178f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
1179f32d5d43SBarry Smith     }
1180f32d5d43SBarry Smith     b1 += bm4;
1181f32d5d43SBarry Smith     b2 += bm4;
1182f32d5d43SBarry Smith     b3 += bm4;
1183f32d5d43SBarry Smith     b4 += bm4;
1184f32d5d43SBarry Smith   }
1185f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
1186f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1187f32d5d43SBarry Smith       r1 = 0.0;
1188f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
1189f32d5d43SBarry Smith       aj  = a->j + a->i[i];
1190f32d5d43SBarry Smith       aa  = a->a + a->i[i];
1191f32d5d43SBarry Smith 
1192f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1193f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
1194f32d5d43SBarry Smith       }
1195f32d5d43SBarry Smith       c[col*am + i]     = r1;
1196f32d5d43SBarry Smith     }
1197f32d5d43SBarry Smith     b1 += bm;
1198f32d5d43SBarry Smith   }
1199dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
12008c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
12018c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1202f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1203f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1204f32d5d43SBarry Smith   PetscFunctionReturn(0);
1205f32d5d43SBarry Smith }
1206f32d5d43SBarry Smith 
1207f32d5d43SBarry Smith /*
12084324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1209f32d5d43SBarry Smith */
1210f32d5d43SBarry Smith #undef __FUNCT__
1211f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
1212f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1213f32d5d43SBarry Smith {
1214f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1215f32d5d43SBarry Smith   PetscErrorCode ierr;
1216dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1217dd6ea824SBarry Smith   MatScalar      *aa;
1218d0f46423SBarry 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;
12194324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1220f32d5d43SBarry Smith 
1221f32d5d43SBarry Smith   PetscFunctionBegin;
1222f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
12238c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
12248c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1225f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12264324174eSBarry Smith 
12274324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
12284324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
12294324174eSBarry Smith       colam = col*am;
12304324174eSBarry Smith       arm   = a->compressedrow.nrows;
12314324174eSBarry Smith       ii    = a->compressedrow.i;
12324324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12334324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12344324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
12354324174eSBarry Smith 	n   = ii[i+1] - ii[i];
12364324174eSBarry Smith 	aj  = a->j + ii[i];
12374324174eSBarry Smith 	aa  = a->a + ii[i];
12384324174eSBarry Smith 	for (j=0; j<n; j++) {
12394324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
12404324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
12414324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
12424324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
12434324174eSBarry Smith 	}
12444324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
12454324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
12464324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
12474324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
12484324174eSBarry Smith       }
12494324174eSBarry Smith       b1 += bm4;
12504324174eSBarry Smith       b2 += bm4;
12514324174eSBarry Smith       b3 += bm4;
12524324174eSBarry Smith       b4 += bm4;
12534324174eSBarry Smith     }
12544324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
12554324174eSBarry Smith       colam = col*am;
12564324174eSBarry Smith       arm   = a->compressedrow.nrows;
12574324174eSBarry Smith       ii    = a->compressedrow.i;
12584324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12594324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12604324174eSBarry Smith 	r1 = 0.0;
12614324174eSBarry Smith 	n   = ii[i+1] - ii[i];
12624324174eSBarry Smith 	aj  = a->j + ii[i];
12634324174eSBarry Smith 	aa  = a->a + ii[i];
12644324174eSBarry Smith 
12654324174eSBarry Smith 	for (j=0; j<n; j++) {
12664324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
12674324174eSBarry Smith 	}
1268a2ea699eSBarry Smith 	c[colam + ridx[i]] += r1;
12694324174eSBarry Smith       }
12704324174eSBarry Smith       b1 += bm;
12714324174eSBarry Smith     }
12724324174eSBarry Smith   } else {
1273f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
1274f32d5d43SBarry Smith       colam = col*am;
1275f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1276f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
1277f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1278f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1279f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1280f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1281f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
1282f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
1283f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
1284f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
1285f32d5d43SBarry Smith 	}
1286f32d5d43SBarry Smith 	c[colam + i]       += r1;
1287f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
1288f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
1289f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
1290f32d5d43SBarry Smith       }
1291f32d5d43SBarry Smith       b1 += bm4;
1292f32d5d43SBarry Smith       b2 += bm4;
1293f32d5d43SBarry Smith       b3 += bm4;
1294f32d5d43SBarry Smith       b4 += bm4;
1295f32d5d43SBarry Smith     }
1296f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
1297a2ea699eSBarry Smith       colam = col*am;
1298f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1299f32d5d43SBarry Smith 	r1 = 0.0;
1300f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
1301f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
1302f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
1303f32d5d43SBarry Smith 
1304f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
1305f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
1306f32d5d43SBarry Smith 	}
1307a2ea699eSBarry Smith 	c[colam + i]     += r1;
1308f32d5d43SBarry Smith       }
1309f32d5d43SBarry Smith       b1 += bm;
1310f32d5d43SBarry Smith     }
13114324174eSBarry Smith   }
1312dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13138c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
13148c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13159123193aSHong Zhang   PetscFunctionReturn(0);
13169123193aSHong Zhang }
1317b1683b59SHong Zhang 
1318b1683b59SHong Zhang #undef __FUNCT__
1319b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
1320b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1321c8db22f6SHong Zhang {
1322c8db22f6SHong Zhang   PetscErrorCode ierr;
13232f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
13242f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13252f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
13262f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13272f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13282f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1329c8db22f6SHong Zhang 
1330c8db22f6SHong Zhang   PetscFunctionBegin;
13312f5f1f90SHong Zhang   btval_den=btdense->v;
13322f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13332f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13342f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13352f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1336d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13372f5f1f90SHong Zhang       btcol = bj + bi[col];
13382f5f1f90SHong Zhang       btval = ba + bi[col];
13392f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1340d2853540SHong Zhang       for (j=0; j<anz; j++){
13412f5f1f90SHong Zhang         brow            = btcol[j];
13422f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1343c8db22f6SHong Zhang       }
1344c8db22f6SHong Zhang     }
13452f5f1f90SHong Zhang     btval_den += m;
1346c8db22f6SHong Zhang   }
1347c8db22f6SHong Zhang   PetscFunctionReturn(0);
1348c8db22f6SHong Zhang }
1349c8db22f6SHong Zhang 
1350c8db22f6SHong Zhang #undef __FUNCT__
1351b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
1352b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13538972f759SHong Zhang {
13548972f759SHong Zhang   PetscErrorCode ierr;
1355b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
13562f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
1357b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
13582f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
13598972f759SHong Zhang 
13608972f759SHong Zhang   PetscFunctionBegin;
13618972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
1362a3fe58edSHong Zhang   ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1363b2d2b04fSHong Zhang   cp_den = ca_den;
13642f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13652f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
13662f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
13672f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
13682f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
13692f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
1370b2d2b04fSHong Zhang     }
1371b2d2b04fSHong Zhang     cp_den += m;
1372b2d2b04fSHong Zhang   }
1373a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
13748972f759SHong Zhang   PetscFunctionReturn(0);
13758972f759SHong Zhang }
13768972f759SHong Zhang 
1377e99be685SHong Zhang /*
1378e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
1379e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
1380e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
1381e99be685SHong Zhang  */
1382e99be685SHong Zhang #undef __FUNCT__
1383e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
13841a83f524SJed Brown PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1385e99be685SHong Zhang {
1386e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1387e99be685SHong Zhang   PetscErrorCode ierr;
1388e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
1389e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
1390e99be685SHong Zhang   PetscInt       *cspidx;
1391e99be685SHong Zhang 
1392e99be685SHong Zhang   PetscFunctionBegin;
1393e99be685SHong Zhang   *nn = n;
1394e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1395e99be685SHong Zhang   if (symmetric) {
1396e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
13971a83f524SJed Brown     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr);
1398e99be685SHong Zhang   } else {
1399e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
1400e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1401e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
1402e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
1403e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
1404e99be685SHong Zhang     jj = a->j;
1405e99be685SHong Zhang     for (i=0; i<nz; i++) {
1406e99be685SHong Zhang       collengths[jj[i]]++;
1407e99be685SHong Zhang     }
1408e99be685SHong Zhang     cia[0] = oshift;
1409e99be685SHong Zhang     for (i=0; i<n; i++) {
1410e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
1411e99be685SHong Zhang     }
1412e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
1413e99be685SHong Zhang     jj   = a->j;
1414e99be685SHong Zhang     for (row=0; row<m; row++) {
1415e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
1416e99be685SHong Zhang       for (i=0; i<mr; i++) {
1417e99be685SHong Zhang         col = *jj++;
1418e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1419e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1420e99be685SHong Zhang       }
1421e99be685SHong Zhang     }
1422e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
1423e99be685SHong Zhang     *ia = cia; *ja = cja;
1424e99be685SHong Zhang     *spidx = cspidx;
1425e99be685SHong Zhang   }
1426e99be685SHong Zhang   PetscFunctionReturn(0);
1427e99be685SHong Zhang }
1428e99be685SHong Zhang 
1429e99be685SHong Zhang #undef __FUNCT__
1430e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
14311a83f524SJed Brown PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1432e99be685SHong Zhang {
1433e99be685SHong Zhang   PetscErrorCode ierr;
1434e99be685SHong Zhang 
1435e99be685SHong Zhang   PetscFunctionBegin;
1436e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
1437e99be685SHong Zhang 
1438e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
1439e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
1440e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1441e99be685SHong Zhang   PetscFunctionReturn(0);
1442e99be685SHong Zhang }
1443e99be685SHong Zhang 
14448972f759SHong Zhang #undef __FUNCT__
1445b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1446b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1447b1683b59SHong Zhang {
1448b1683b59SHong Zhang   PetscErrorCode ierr;
14491a83f524SJed Brown   PetscInt       i,n,nrows,N,j,k,m,ncols,col,cm;
14501a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1451b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1452b1683b59SHong Zhang   IS             *isa;
1453b1683b59SHong Zhang   PetscBool      flg1,flg2;
1454b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1455e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1456d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1457b1683b59SHong Zhang 
1458b1683b59SHong Zhang   PetscFunctionBegin;
1459b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1460e99be685SHong Zhang 
1461b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1462251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1463251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1464b1683b59SHong Zhang   if (flg1 || flg2) {
1465b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1466b1683b59SHong Zhang   }
1467b1683b59SHong Zhang 
1468b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1469b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1470b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1471b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1472b1683b59SHong Zhang   c->rstart = 0;
1473b1683b59SHong Zhang 
1474b1683b59SHong Zhang   c->ncolors = nis;
1475b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1476b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1477d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1478d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1479d2853540SHong Zhang   colorforrow[0]    = 0;
1480d2853540SHong Zhang   rows_i            = rows;
1481d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1482b1683b59SHong Zhang 
1483d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1484d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1485d2853540SHong Zhang   colorforcol[0] = 0;
1486d2853540SHong Zhang   columns_i      = columns;
1487d2853540SHong Zhang 
1488f2691f8dSHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); /* column-wise storage of mat */
1489b1683b59SHong Zhang 
1490b28d80bdSHong Zhang   cm = c->m;
1491b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1492e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1493b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1494b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1495b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1496b1683b59SHong Zhang     c->ncolumns[i] = n;
1497b1683b59SHong Zhang     if (n) {
1498d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1499b1683b59SHong Zhang     }
1500d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1501d2853540SHong Zhang     columns_i       += n;
1502b1683b59SHong Zhang 
1503b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1504e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1505e99be685SHong Zhang 
1506b1683b59SHong Zhang     /* loop over columns*/
1507b1683b59SHong Zhang     for (j=0; j<n; j++) {
1508b1683b59SHong Zhang       col     = is[j];
1509d2853540SHong Zhang       row_idx = cj + ci[col];
1510b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1511b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1512b1683b59SHong Zhang       for (k=0; k<m; k++) {
1513e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1514d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1515b1683b59SHong Zhang       }
1516b1683b59SHong Zhang     }
1517b1683b59SHong Zhang     /* count the number of hits */
1518b1683b59SHong Zhang     nrows = 0;
1519e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1520b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1521b1683b59SHong Zhang     }
1522b1683b59SHong Zhang     c->nrows[i]      = nrows;
1523d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1524d2853540SHong Zhang 
1525b1683b59SHong Zhang     nrows       = 0;
1526e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1527b1683b59SHong Zhang       if (rowhit[j]) {
1528d2853540SHong Zhang         rows_i[nrows]            = j;
1529e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1530b1683b59SHong Zhang         nrows++;
1531b1683b59SHong Zhang       }
1532b1683b59SHong Zhang     }
1533b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1534d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1535b1683b59SHong Zhang   }
1536f2691f8dSHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr);
1537b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1538b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1539d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1540d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1541d2853540SHong Zhang #endif
1542b28d80bdSHong Zhang 
1543d2853540SHong Zhang   c->colorforrow     = colorforrow;
1544d2853540SHong Zhang   c->rows            = rows;
1545d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1546d2853540SHong Zhang   c->colorforcol     = colorforcol;
1547d2853540SHong Zhang   c->columns         = columns;
1548f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1549b1683b59SHong Zhang   PetscFunctionReturn(0);
1550b1683b59SHong Zhang }
1551