xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision b561aa9d943d3781294c5a9d05dba64bef806d4b)
1be1d678aSKris Buschelman 
2d50806bdSBarry Smith /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4d50806bdSBarry Smith           C = A * B
5d50806bdSBarry Smith */
6d50806bdSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
116284ec50SHong Zhang 
12ec01deb9SMatthew Knepley EXTERN_C_BEGIN
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;
185c66b693SKris Buschelman 
195c66b693SKris Buschelman   PetscFunctionBegin;
2026be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
2126be0446SHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
2226be0446SHong Zhang   }
2326be0446SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
245c66b693SKris Buschelman   PetscFunctionReturn(0);
255c66b693SKris Buschelman }
26ec01deb9SMatthew Knepley EXTERN_C_END
271c24bd37SHong Zhang 
2826be0446SHong Zhang #undef __FUNCT__
29b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ"
30b561aa9dSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
3158c24d83SHong Zhang {
32dfbe8321SBarry Smith   PetscErrorCode     ierr;
33a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
34b561aa9dSHong Zhang   PetscInt           *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj;
35b561aa9dSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0;
36be0fcf8dSHong Zhang   PetscBT            lnkbt;
3758c24d83SHong Zhang 
3858c24d83SHong Zhang   PetscFunctionBegin;
3958c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
4058c24d83SHong Zhang   /* free space for accumulating nonzero column info */
4138baddfdSBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4258c24d83SHong Zhang   ci[0] = 0;
4358c24d83SHong Zhang 
44be0fcf8dSHong Zhang   /* create and initialize a linked list */
45be0fcf8dSHong Zhang   nlnk = bn+1;
46be0fcf8dSHong Zhang   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
4758c24d83SHong Zhang 
48c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
49a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5058c24d83SHong Zhang   current_space = free_space;
5158c24d83SHong Zhang 
5258c24d83SHong Zhang   /* Determine symbolic info for each row of the product: */
5358c24d83SHong Zhang   for (i=0;i<am;i++) {
5458c24d83SHong Zhang     anzi = ai[i+1] - ai[i];
5558c24d83SHong Zhang     cnzi = 0;
562d09714cSHong Zhang     j    = anzi;
57b561aa9dSHong Zhang     aj   = Aj + ai[i];
582d09714cSHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
592d09714cSHong Zhang       j--;
60b561aa9dSHong Zhang       brow = aj[j];
6158c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
6258c24d83SHong Zhang       bjj  = bj + bi[brow];
631c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
64be0fcf8dSHong Zhang       ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
651c239cc6SHong Zhang       cnzi += nlnk;
6658c24d83SHong Zhang     }
6758c24d83SHong Zhang 
6858c24d83SHong Zhang     /* If free space is not available, make more free space */
6958c24d83SHong Zhang     /* Double the amount of total space in the list */
7058c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
714238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
72b561aa9dSHong Zhang       ndouble++;
7358c24d83SHong Zhang     }
7458c24d83SHong Zhang 
75c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
76be0fcf8dSHong Zhang     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
77c5db241fSHong Zhang     current_space->array           += cnzi;
7858c24d83SHong Zhang     current_space->local_used      += cnzi;
7958c24d83SHong Zhang     current_space->local_remaining -= cnzi;
8058c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
8158c24d83SHong Zhang   }
8258c24d83SHong Zhang 
8358c24d83SHong Zhang   /* Column indices are in the list of free space */
8458c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
8558c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
8638baddfdSBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
87a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
88be0fcf8dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
8958c24d83SHong Zhang 
90b561aa9dSHong Zhang   *Ci           = ci;
91b561aa9dSHong Zhang   *Cj           = cj;
92b561aa9dSHong Zhang   *nspacedouble = ndouble;
93b561aa9dSHong Zhang   PetscFunctionReturn(0);
94b561aa9dSHong Zhang }
95b561aa9dSHong Zhang 
96b561aa9dSHong Zhang #undef __FUNCT__
97b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
98b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
99b561aa9dSHong Zhang {
100b561aa9dSHong Zhang   PetscErrorCode ierr;
101b561aa9dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
102b561aa9dSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj;
103b561aa9dSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
104b561aa9dSHong Zhang   MatScalar      *ca;
105b561aa9dSHong Zhang   PetscReal      afill;
106b561aa9dSHong Zhang 
107b561aa9dSHong Zhang   PetscFunctionBegin;
108b561aa9dSHong Zhang   /* Get ci and cj */
109b561aa9dSHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
110b561aa9dSHong Zhang 
11158c24d83SHong Zhang   /* Allocate space for ca */
11258c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
11358c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
11458c24d83SHong Zhang 
11526be0446SHong Zhang   /* put together the new symbolic matrix */
1167adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
11758c24d83SHong Zhang 
11858c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
11958c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
12058c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
121e6b907acSBarry Smith   c->free_a   = PETSC_TRUE;
122e6b907acSBarry Smith   c->free_ij  = PETSC_TRUE;
12358c24d83SHong Zhang   c->nonew    = 0;
12458c24d83SHong Zhang 
1257212da7cSHong Zhang   /* set MatInfo */
1267212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
127f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1287212da7cSHong Zhang   c->maxnz                     = ci[am];
1297212da7cSHong Zhang   c->nz                        = ci[am];
1307212da7cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
1317212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1327212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1337212da7cSHong Zhang 
1347212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1357212da7cSHong Zhang   if (ci[am]) {
136f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
137f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
138f2b054eeSHong Zhang   } else {
139f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
140be0fcf8dSHong Zhang   }
141f2b054eeSHong Zhang #endif
14258c24d83SHong Zhang   PetscFunctionReturn(0);
14358c24d83SHong Zhang }
144d50806bdSBarry Smith 
14526be0446SHong Zhang #undef __FUNCT__
14626be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
147dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
148d50806bdSBarry Smith {
149dfbe8321SBarry Smith   PetscErrorCode ierr;
150d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
151d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
152d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
153d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
15438baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
155d0f46423SBarry Smith   PetscInt       am=A->rmap->N,cm=C->rmap->N;
156c124e916SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow,nextb;
157c124e916SHong Zhang   MatScalar      *aa=a->a,*ba=b->a,*baj,*ca=c->a;
158d50806bdSBarry Smith 
159d50806bdSBarry Smith   PetscFunctionBegin;
160c124e916SHong Zhang   /* clean old values in C */
161c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
162d50806bdSBarry Smith   /* Traverse A row-wise. */
163d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
164d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
165d50806bdSBarry Smith   for (i=0;i<am;i++) {
166d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
167d50806bdSBarry Smith     for (j=0;j<anzi;j++) {
168d50806bdSBarry Smith       brow = *aj++;
169d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
170d50806bdSBarry Smith       bjj  = bj + bi[brow];
171d50806bdSBarry Smith       baj  = ba + bi[brow];
172c124e916SHong Zhang       nextb = 0;
173c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
174c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
175c124e916SHong Zhang           ca[k] += (*aa)*baj[nextb++];
176c124e916SHong Zhang         }
177d50806bdSBarry Smith       }
178d50806bdSBarry Smith       flops += 2*bnzi;
179d50806bdSBarry Smith       aa++;
180d50806bdSBarry Smith     }
181d50806bdSBarry Smith     cnzi = ci[i+1] - ci[i];
182d50806bdSBarry Smith     ca += cnzi;
183d50806bdSBarry Smith     cj += cnzi;
184d50806bdSBarry Smith   }
185716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
186716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
187716bacf3SKris Buschelman 
188d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
189d50806bdSBarry Smith   PetscFunctionReturn(0);
190d50806bdSBarry Smith }
191bc011b1eSHong Zhang 
192d2853540SHong Zhang /* This routine is not used. Should be removed! */
193bc011b1eSHong Zhang #undef __FUNCT__
1946fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
1956fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1965df89d91SHong Zhang {
197bc011b1eSHong Zhang   PetscErrorCode ierr;
198bc011b1eSHong Zhang 
199bc011b1eSHong Zhang   PetscFunctionBegin;
200bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
2016fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
202bc011b1eSHong Zhang   }
2036fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
204bc011b1eSHong Zhang   PetscFunctionReturn(0);
205bc011b1eSHong Zhang }
206bc011b1eSHong Zhang 
207bc011b1eSHong Zhang #undef __FUNCT__
208e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
209e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
2102128a86cSHong Zhang {
2112128a86cSHong Zhang   PetscErrorCode      ierr;
212e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
2132128a86cSHong Zhang 
2142128a86cSHong Zhang   PetscFunctionBegin;
215b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
2162128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
2172128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
2182128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
2192128a86cSHong Zhang   PetscFunctionReturn(0);
2202128a86cSHong Zhang }
2212128a86cSHong Zhang 
2222128a86cSHong Zhang #undef __FUNCT__
2232128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
2242128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
2252128a86cSHong Zhang {
2262128a86cSHong Zhang   PetscErrorCode      ierr;
2272128a86cSHong Zhang   PetscContainer      container;
228e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
2292128a86cSHong Zhang 
2302128a86cSHong Zhang   PetscFunctionBegin;
231e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
2322128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
2332128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
2342128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
2352128a86cSHong Zhang   if (A->ops->destroy) {
2362128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2372128a86cSHong Zhang   }
238e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
2392128a86cSHong Zhang   PetscFunctionReturn(0);
2402128a86cSHong Zhang }
2412128a86cSHong Zhang 
2422128a86cSHong Zhang #undef __FUNCT__
2436fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
2446fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
245bc011b1eSHong Zhang {
2465df89d91SHong Zhang   PetscErrorCode      ierr;
24781d82fe4SHong Zhang   Mat                 Bt;
24881d82fe4SHong Zhang   PetscInt            *bti,*btj;
249e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
2502128a86cSHong Zhang   PetscContainer      container;
251d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
252d2853540SHong Zhang 
25381d82fe4SHong Zhang   PetscFunctionBegin;
254d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
25581d82fe4SHong Zhang    /* create symbolic Bt */
25681d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
25781d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
25881d82fe4SHong Zhang 
25981d82fe4SHong Zhang   /* get symbolic C=A*Bt */
26081d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
26181d82fe4SHong Zhang 
2622128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
263e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
2642128a86cSHong Zhang 
2652128a86cSHong Zhang   /* attach the supporting struct to C */
2662128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
2672128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
268e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
269e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
2702128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
2712128a86cSHong Zhang 
2722128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
2732128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
2742128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
2752128a86cSHong Zhang 
276d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
277d2853540SHong Zhang   etime2 += tf - t0;
278d2853540SHong Zhang 
279f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
2802128a86cSHong Zhang   if (multtrans->usecoloring){
281b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
282b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
2832128a86cSHong Zhang     ISColoring           iscoloring;
2842128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
285d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
286e8354b3eSHong Zhang 
287e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
288502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
289d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
290d2853540SHong Zhang     etime0 += tf - t0;
291d2853540SHong Zhang 
292d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
293b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
2942128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
2952128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
296e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
297d2853540SHong Zhang     etime01 += tf - t0;
2982128a86cSHong Zhang 
299e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
3002128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
3012128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
3022128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
3032128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
3042128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
3052128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
3062128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
3072128a86cSHong Zhang 
3082128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
3092128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
3102128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
3112128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
3122128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
3132128a86cSHong Zhang     multtrans->ABt_den = C_dense;
314e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
315e8354b3eSHong Zhang     etime1 += tf - t0;
316f94ccd6cSHong Zhang 
317f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
318f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
319f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));
320f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
321f94ccd6cSHong Zhang #endif
3222128a86cSHong Zhang   }
323*e99be685SHong Zhang   /* clean up */
324*e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
325*e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
3262128a86cSHong Zhang 
327f94ccd6cSHong Zhang 
328f94ccd6cSHong Zhang 
32981d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
33081d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
3311fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
3321fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
3331fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
3341fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
3351fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
3361fa1209cSHong Zhang   MatScalar          *ca;
3371fa1209cSHong Zhang   PetscBT            lnkbt;
3381fa1209cSHong Zhang   PetscReal          afill;
3395df89d91SHong Zhang 
3401fa1209cSHong Zhang   /* Allocate row pointer array ci  */
3411fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
3421fa1209cSHong Zhang   ci[0] = 0;
3431fa1209cSHong Zhang 
3441fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
3451fa1209cSHong Zhang   nlnk = bm+1;
3461fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3471fa1209cSHong Zhang 
3481fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
3491fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
3501fa1209cSHong Zhang   current_space = free_space;
3511fa1209cSHong Zhang 
3521fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
3531fa1209cSHong Zhang   for (i=0; i<am; i++) {
3541fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
3551fa1209cSHong Zhang     cnzi = 0;
3561fa1209cSHong Zhang     acol = aj + ai[i];
3571fa1209cSHong Zhang     for (j=0; j<bm; j++){
3581fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
3591fa1209cSHong Zhang       bcol= bj + bi[j];
36081d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
3611fa1209cSHong Zhang       ka = 0; kb = 0;
3621fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
36381d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
36481d82fe4SHong Zhang         if (ka == anzi) break;
36581d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
36681d82fe4SHong Zhang         if (kb == bnzj) break;
3671fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
3681fa1209cSHong Zhang           index[0] = j;
3691fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3701fa1209cSHong Zhang           cnzi++;
3711fa1209cSHong Zhang           break;
3721fa1209cSHong Zhang         }
3731fa1209cSHong Zhang       }
3741fa1209cSHong Zhang     }
3751fa1209cSHong Zhang 
3761fa1209cSHong Zhang     /* If free space is not available, make more free space */
3771fa1209cSHong Zhang     /* Double the amount of total space in the list */
3781fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
3791fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
3801fa1209cSHong Zhang       nspacedouble++;
3811fa1209cSHong Zhang     }
3821fa1209cSHong Zhang 
3831fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
3841fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3851fa1209cSHong Zhang     current_space->array           += cnzi;
3861fa1209cSHong Zhang     current_space->local_used      += cnzi;
3871fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
3881fa1209cSHong Zhang 
3891fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
3901fa1209cSHong Zhang   }
3911fa1209cSHong Zhang 
3921fa1209cSHong Zhang 
3931fa1209cSHong Zhang   /* Column indices are in the list of free space.
3941fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
3951fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
3961fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3971fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3981fa1209cSHong Zhang 
3991fa1209cSHong Zhang   /* Allocate space for ca */
4001fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4011fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4021fa1209cSHong Zhang 
4031fa1209cSHong Zhang   /* put together the new symbolic matrix */
4041fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
4051fa1209cSHong Zhang 
4061fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4071fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
4081fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
4091fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
4101fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
4111fa1209cSHong Zhang   c->nonew    = 0;
4121fa1209cSHong Zhang 
4131fa1209cSHong Zhang   /* set MatInfo */
4141fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
4151fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
4161fa1209cSHong Zhang   c->maxnz                     = ci[am];
4171fa1209cSHong Zhang   c->nz                        = ci[am];
4181fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
4191fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
4201fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
4211fa1209cSHong Zhang 
4221fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
4231fa1209cSHong Zhang   if (ci[am]) {
4241fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
4256fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
4261fa1209cSHong Zhang   } else {
4271fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
4281fa1209cSHong Zhang   }
4291fa1209cSHong Zhang #endif
43081d82fe4SHong Zhang #endif
4315df89d91SHong Zhang   PetscFunctionReturn(0);
4325df89d91SHong Zhang }
4335df89d91SHong Zhang 
4346973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
4355df89d91SHong Zhang #undef __FUNCT__
4366fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
4376fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
4385df89d91SHong Zhang {
4395df89d91SHong Zhang   PetscErrorCode ierr;
4405df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
441e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
44281d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
4435df89d91SHong Zhang   PetscLogDouble flops=0.0;
4446973769bSHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
445e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
4462128a86cSHong Zhang   PetscContainer      container;
4476973769bSHong Zhang #if defined(USE_ARRAY)
4486973769bSHong Zhang   MatScalar      *spdot;
4496973769bSHong Zhang #endif
4505df89d91SHong Zhang 
4515df89d91SHong Zhang   PetscFunctionBegin;
452e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
4532128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
4542128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
4552128a86cSHong Zhang   if (multtrans->usecoloring){
456b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
457c8db22f6SHong Zhang     Mat                   Bt_dense;
458c8db22f6SHong Zhang     PetscInt              m,n;
459b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
460a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
461a0b95ffbSSatish Balay 
4622128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
463c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
464c8db22f6SHong Zhang 
465b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
4662128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
467b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
4682128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
469b2d2b04fSHong Zhang     etime0 += tf - t0;
470c8db22f6SHong Zhang 
471c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
4722128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
4732128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
4742128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
4752128a86cSHong Zhang     etime2 += tf - t0;
476c8db22f6SHong Zhang 
4772128a86cSHong Zhang     /* Recover C from C_dense */
4782128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
479b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
4802128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
4812128a86cSHong Zhang     etime1 += tf - t0;
482f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
483f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
484f94ccd6cSHong Zhang #endif
485c8db22f6SHong Zhang     PetscFunctionReturn(0);
486c8db22f6SHong Zhang   }
487c8db22f6SHong Zhang 
4886973769bSHong Zhang #if defined(USE_ARRAY)
4896973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
490e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
491e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
4926973769bSHong Zhang #endif
4936973769bSHong Zhang 
49481d82fe4SHong Zhang   /* clear old values in C */
49581d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
4965df89d91SHong Zhang 
4971fa1209cSHong Zhang   for (i=0; i<cm; i++) {
49881d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
49981d82fe4SHong Zhang     acol = aj + ai[i];
5006973769bSHong Zhang     aval = aa + ai[i];
5011fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
5021fa1209cSHong Zhang     ccol = cj + ci[i];
5036973769bSHong Zhang     cval = ca + ci[i];
5041fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
50581d82fe4SHong Zhang       brow = ccol[j];
50681d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
50781d82fe4SHong Zhang       bcol = bj + bi[brow];
5086973769bSHong Zhang       bval = ba + bi[brow];
5096973769bSHong Zhang 
51081d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
5116973769bSHong Zhang #if defined(USE_ARRAY)
5126973769bSHong Zhang       /* put ba to spdot array */
5136973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
5146973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
5156973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
5166973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
5176973769bSHong Zhang       }
5186973769bSHong Zhang       /* zero spdot array */
5196973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
5206973769bSHong Zhang #else
52181d82fe4SHong Zhang       nexta = 0; nextb = 0;
52281d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
52381d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
52481d82fe4SHong Zhang         if (nexta == anzi) break;
52581d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
52681d82fe4SHong Zhang         if (nextb == bnzj) break;
52781d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
5286973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
52981d82fe4SHong Zhang           nexta++; nextb++;
53081d82fe4SHong Zhang           flops += 2;
5311fa1209cSHong Zhang         }
5321fa1209cSHong Zhang       }
5336973769bSHong Zhang #endif
53481d82fe4SHong Zhang     }
53581d82fe4SHong Zhang   }
53681d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
53781d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
53881d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
5396973769bSHong Zhang #if defined(USE_ARRAY)
5406973769bSHong Zhang   ierr = PetscFree(spdot);
5416973769bSHong Zhang #endif
5425df89d91SHong Zhang   PetscFunctionReturn(0);
5435df89d91SHong Zhang }
5445df89d91SHong Zhang 
5455df89d91SHong Zhang #undef __FUNCT__
54675648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
54775648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
5485df89d91SHong Zhang   PetscErrorCode ierr;
5495df89d91SHong Zhang 
5505df89d91SHong Zhang   PetscFunctionBegin;
5515df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
55275648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
5535df89d91SHong Zhang   }
55475648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
5555df89d91SHong Zhang   PetscFunctionReturn(0);
5565df89d91SHong Zhang }
5575df89d91SHong Zhang 
5585df89d91SHong Zhang #undef __FUNCT__
55975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
56075648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
5615df89d91SHong Zhang {
562bc011b1eSHong Zhang   PetscErrorCode ierr;
563bc011b1eSHong Zhang   Mat            At;
56438baddfdSBarry Smith   PetscInt       *ati,*atj;
565bc011b1eSHong Zhang 
566bc011b1eSHong Zhang   PetscFunctionBegin;
567bc011b1eSHong Zhang   /* create symbolic At */
568bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
569d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
570bc011b1eSHong Zhang 
571bc011b1eSHong Zhang   /* get symbolic C=At*B */
572bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
573bc011b1eSHong Zhang 
574bc011b1eSHong Zhang   /* clean up */
5756bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
576bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
577bc011b1eSHong Zhang   PetscFunctionReturn(0);
578bc011b1eSHong Zhang }
579bc011b1eSHong Zhang 
580bc011b1eSHong Zhang #undef __FUNCT__
58175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
58275648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
583bc011b1eSHong Zhang {
584bc011b1eSHong Zhang   PetscErrorCode ierr;
5850fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
586d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
587d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
588d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
5890fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
590bc011b1eSHong Zhang 
591bc011b1eSHong Zhang   PetscFunctionBegin;
592bc011b1eSHong Zhang   /* clear old values in C */
593bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
594bc011b1eSHong Zhang 
595bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
596bc011b1eSHong Zhang   for (i=0;i<am;i++) {
597bc011b1eSHong Zhang     bj   = b->j + bi[i];
598bc011b1eSHong Zhang     ba   = b->a + bi[i];
599bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
600bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
601bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
602bc011b1eSHong Zhang       nextb = 0;
6030fbc74f4SHong Zhang       crow  = *aj++;
604bc011b1eSHong Zhang       cjj   = cj + ci[crow];
605bc011b1eSHong Zhang       caj   = ca + ci[crow];
606bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
607bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
6080fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
6090fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
610bc011b1eSHong Zhang           nextb++;
611bc011b1eSHong Zhang         }
612bc011b1eSHong Zhang       }
613bc011b1eSHong Zhang       flops += 2*bnzi;
6140fbc74f4SHong Zhang       aa++;
615bc011b1eSHong Zhang     }
616bc011b1eSHong Zhang   }
617bc011b1eSHong Zhang 
618bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
619bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
620bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
621bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
622bc011b1eSHong Zhang   PetscFunctionReturn(0);
623bc011b1eSHong Zhang }
6249123193aSHong Zhang 
625ec01deb9SMatthew Knepley EXTERN_C_BEGIN
6269123193aSHong Zhang #undef __FUNCT__
6279123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
6289123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
6299123193aSHong Zhang {
6309123193aSHong Zhang   PetscErrorCode ierr;
6319123193aSHong Zhang 
6329123193aSHong Zhang   PetscFunctionBegin;
6339123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6349123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
6359123193aSHong Zhang   }
6369123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
6379123193aSHong Zhang   PetscFunctionReturn(0);
6389123193aSHong Zhang }
639ec01deb9SMatthew Knepley EXTERN_C_END
640edd81eecSMatthew Knepley 
6419123193aSHong Zhang #undef __FUNCT__
6429123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
6439123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
6449123193aSHong Zhang {
6459123193aSHong Zhang   PetscErrorCode ierr;
6469123193aSHong Zhang 
6479123193aSHong Zhang   PetscFunctionBegin;
6485a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
6499123193aSHong Zhang   PetscFunctionReturn(0);
6509123193aSHong Zhang }
6519123193aSHong Zhang 
6529123193aSHong Zhang #undef __FUNCT__
6539123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
6549123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
6559123193aSHong Zhang {
656f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
6579123193aSHong Zhang   PetscErrorCode ierr;
658dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
659dd6ea824SBarry Smith   MatScalar      *aa;
660d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
661f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
6629123193aSHong Zhang 
6639123193aSHong Zhang   PetscFunctionBegin;
664f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
665e32f2f54SBarry 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);
666e32f2f54SBarry 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);
667e32f2f54SBarry 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);
668f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
669f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
670f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
671f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
672f32d5d43SBarry Smith     colam = col*am;
673f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
674f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
675f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
676f32d5d43SBarry Smith       aj  = a->j + a->i[i];
677f32d5d43SBarry Smith       aa  = a->a + a->i[i];
678f32d5d43SBarry Smith       for (j=0; j<n; j++) {
679f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
680f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
681f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
682f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
6839123193aSHong Zhang       }
684f32d5d43SBarry Smith       c[colam + i]       = r1;
685f32d5d43SBarry Smith       c[colam + am + i]  = r2;
686f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
687f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
688f32d5d43SBarry Smith     }
689f32d5d43SBarry Smith     b1 += bm4;
690f32d5d43SBarry Smith     b2 += bm4;
691f32d5d43SBarry Smith     b3 += bm4;
692f32d5d43SBarry Smith     b4 += bm4;
693f32d5d43SBarry Smith   }
694f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
695f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
696f32d5d43SBarry Smith       r1 = 0.0;
697f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
698f32d5d43SBarry Smith       aj  = a->j + a->i[i];
699f32d5d43SBarry Smith       aa  = a->a + a->i[i];
700f32d5d43SBarry Smith 
701f32d5d43SBarry Smith       for (j=0; j<n; j++) {
702f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
703f32d5d43SBarry Smith       }
704f32d5d43SBarry Smith       c[col*am + i]     = r1;
705f32d5d43SBarry Smith     }
706f32d5d43SBarry Smith     b1 += bm;
707f32d5d43SBarry Smith   }
708dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
709f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
710f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
711f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
712f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
713f32d5d43SBarry Smith   PetscFunctionReturn(0);
714f32d5d43SBarry Smith }
715f32d5d43SBarry Smith 
716f32d5d43SBarry Smith /*
7174324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
718f32d5d43SBarry Smith */
719f32d5d43SBarry Smith #undef __FUNCT__
720f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
721f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
722f32d5d43SBarry Smith {
723f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
724f32d5d43SBarry Smith   PetscErrorCode ierr;
725dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
726dd6ea824SBarry Smith   MatScalar      *aa;
727d0f46423SBarry 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;
7284324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
729f32d5d43SBarry Smith 
730f32d5d43SBarry Smith   PetscFunctionBegin;
731f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
732f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
733f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
734f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
7354324174eSBarry Smith 
7364324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
7374324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
7384324174eSBarry Smith       colam = col*am;
7394324174eSBarry Smith       arm   = a->compressedrow.nrows;
7404324174eSBarry Smith       ii    = a->compressedrow.i;
7414324174eSBarry Smith       ridx  = a->compressedrow.rindex;
7424324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
7434324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
7444324174eSBarry Smith 	n   = ii[i+1] - ii[i];
7454324174eSBarry Smith 	aj  = a->j + ii[i];
7464324174eSBarry Smith 	aa  = a->a + ii[i];
7474324174eSBarry Smith 	for (j=0; j<n; j++) {
7484324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
7494324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
7504324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
7514324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
7524324174eSBarry Smith 	}
7534324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
7544324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
7554324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
7564324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
7574324174eSBarry Smith       }
7584324174eSBarry Smith       b1 += bm4;
7594324174eSBarry Smith       b2 += bm4;
7604324174eSBarry Smith       b3 += bm4;
7614324174eSBarry Smith       b4 += bm4;
7624324174eSBarry Smith     }
7634324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
7644324174eSBarry Smith       colam = col*am;
7654324174eSBarry Smith       arm   = a->compressedrow.nrows;
7664324174eSBarry Smith       ii    = a->compressedrow.i;
7674324174eSBarry Smith       ridx  = a->compressedrow.rindex;
7684324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
7694324174eSBarry Smith 	r1 = 0.0;
7704324174eSBarry Smith 	n   = ii[i+1] - ii[i];
7714324174eSBarry Smith 	aj  = a->j + ii[i];
7724324174eSBarry Smith 	aa  = a->a + ii[i];
7734324174eSBarry Smith 
7744324174eSBarry Smith 	for (j=0; j<n; j++) {
7754324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
7764324174eSBarry Smith 	}
7774324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
7784324174eSBarry Smith       }
7794324174eSBarry Smith       b1 += bm;
7804324174eSBarry Smith     }
7814324174eSBarry Smith   } else {
782f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
783f32d5d43SBarry Smith       colam = col*am;
784f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
785f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
786f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
787f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
788f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
789f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
790f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
791f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
792f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
793f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
794f32d5d43SBarry Smith 	}
795f32d5d43SBarry Smith 	c[colam + i]       += r1;
796f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
797f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
798f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
799f32d5d43SBarry Smith       }
800f32d5d43SBarry Smith       b1 += bm4;
801f32d5d43SBarry Smith       b2 += bm4;
802f32d5d43SBarry Smith       b3 += bm4;
803f32d5d43SBarry Smith       b4 += bm4;
804f32d5d43SBarry Smith     }
805f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
806f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
807f32d5d43SBarry Smith 	r1 = 0.0;
808f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
809f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
810f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
811f32d5d43SBarry Smith 
812f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
813f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
814f32d5d43SBarry Smith 	}
815f32d5d43SBarry Smith 	c[col*am + i]     += r1;
816f32d5d43SBarry Smith       }
817f32d5d43SBarry Smith       b1 += bm;
818f32d5d43SBarry Smith     }
8194324174eSBarry Smith   }
820dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
821f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
822f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
8239123193aSHong Zhang   PetscFunctionReturn(0);
8249123193aSHong Zhang }
825b1683b59SHong Zhang 
826b1683b59SHong Zhang #undef __FUNCT__
827b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
828b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
829c8db22f6SHong Zhang {
830c8db22f6SHong Zhang   PetscErrorCode ierr;
8312f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
8322f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
8332f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
8342f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
8352f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
8362f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
837c8db22f6SHong Zhang 
838c8db22f6SHong Zhang   PetscFunctionBegin;
8392f5f1f90SHong Zhang   btval_den=btdense->v;
8402f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
8412f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
8422f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
8432f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
844d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
8452f5f1f90SHong Zhang       btcol = bj + bi[col];
8462f5f1f90SHong Zhang       btval = ba + bi[col];
8472f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
848d2853540SHong Zhang       for (j=0; j<anz; j++){
8492f5f1f90SHong Zhang         brow            = btcol[j];
8502f5f1f90SHong Zhang         btval_den[brow] = btval[j];
851c8db22f6SHong Zhang       }
852c8db22f6SHong Zhang     }
8532f5f1f90SHong Zhang     btval_den += m;
854c8db22f6SHong Zhang   }
855c8db22f6SHong Zhang   PetscFunctionReturn(0);
856c8db22f6SHong Zhang }
857c8db22f6SHong Zhang 
858c8db22f6SHong Zhang #undef __FUNCT__
859b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
860b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
8618972f759SHong Zhang {
8628972f759SHong Zhang   PetscErrorCode ierr;
863b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
8642f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
865b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
8662f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
8678972f759SHong Zhang 
8688972f759SHong Zhang   PetscFunctionBegin;
8698972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
870b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
871b2d2b04fSHong Zhang   cp_den = ca_den;
8722f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
8732f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
8742f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
8752f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
8762f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
8772f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
878b2d2b04fSHong Zhang     }
879b2d2b04fSHong Zhang     cp_den += m;
880b2d2b04fSHong Zhang   }
881b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
8828972f759SHong Zhang   PetscFunctionReturn(0);
8838972f759SHong Zhang }
8848972f759SHong Zhang 
885*e99be685SHong Zhang /*
886*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
887*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
888*e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
889*e99be685SHong Zhang  */
890*e99be685SHong Zhang #undef __FUNCT__
891*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
892*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
893*e99be685SHong Zhang {
894*e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
895*e99be685SHong Zhang   PetscErrorCode ierr;
896*e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
897*e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
898*e99be685SHong Zhang   PetscInt       *cspidx;
899*e99be685SHong Zhang 
900*e99be685SHong Zhang   PetscFunctionBegin;
901*e99be685SHong Zhang   *nn = n;
902*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
903*e99be685SHong Zhang   if (symmetric) {
904*e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
905*e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
906*e99be685SHong Zhang   } else {
907*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
908*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
909*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
910*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
911*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
912*e99be685SHong Zhang     jj = a->j;
913*e99be685SHong Zhang     for (i=0; i<nz; i++) {
914*e99be685SHong Zhang       collengths[jj[i]]++;
915*e99be685SHong Zhang     }
916*e99be685SHong Zhang     cia[0] = oshift;
917*e99be685SHong Zhang     for (i=0; i<n; i++) {
918*e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
919*e99be685SHong Zhang     }
920*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
921*e99be685SHong Zhang     jj   = a->j;
922*e99be685SHong Zhang     for (row=0; row<m; row++) {
923*e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
924*e99be685SHong Zhang       for (i=0; i<mr; i++) {
925*e99be685SHong Zhang         col = *jj++;
926*e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
927*e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
928*e99be685SHong Zhang       }
929*e99be685SHong Zhang     }
930*e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
931*e99be685SHong Zhang     *ia = cia; *ja = cja;
932*e99be685SHong Zhang     *spidx = cspidx;
933*e99be685SHong Zhang   }
934*e99be685SHong Zhang   PetscFunctionReturn(0);
935*e99be685SHong Zhang }
936*e99be685SHong Zhang 
937*e99be685SHong Zhang #undef __FUNCT__
938*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
939*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
940*e99be685SHong Zhang {
941*e99be685SHong Zhang   PetscErrorCode ierr;
942*e99be685SHong Zhang 
943*e99be685SHong Zhang   PetscFunctionBegin;
944*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
945*e99be685SHong Zhang 
946*e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
947*e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
948*e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
949*e99be685SHong Zhang   PetscFunctionReturn(0);
950*e99be685SHong Zhang }
951*e99be685SHong Zhang 
9528972f759SHong Zhang #undef __FUNCT__
953b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
954b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
955b1683b59SHong Zhang {
956b1683b59SHong Zhang   PetscErrorCode ierr;
957d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
958b1683b59SHong Zhang   const PetscInt *is;
959b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
960b1683b59SHong Zhang   IS             *isa;
961b1683b59SHong Zhang   PetscBool      done;
962b1683b59SHong Zhang   PetscBool      flg1,flg2;
963b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
964*e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
965d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
966b1683b59SHong Zhang 
967b1683b59SHong Zhang   PetscFunctionBegin;
968b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
969*e99be685SHong Zhang 
970b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
971b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
972b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
973b1683b59SHong Zhang   if (flg1 || flg2) {
974b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
975b1683b59SHong Zhang   }
976b1683b59SHong Zhang 
977b1683b59SHong Zhang   N         = mat->cmap->N/bs;
978b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
979b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
980b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
981b1683b59SHong Zhang   c->rstart = 0;
982b1683b59SHong Zhang 
983b1683b59SHong Zhang   c->ncolors = nis;
984b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
985b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
986d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
987d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
988d2853540SHong Zhang   colorforrow[0]    = 0;
989d2853540SHong Zhang   rows_i            = rows;
990d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
991b1683b59SHong Zhang 
992d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
993d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
994d2853540SHong Zhang   colorforcol[0] = 0;
995d2853540SHong Zhang   columns_i      = columns;
996d2853540SHong Zhang 
997*e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
998b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
999b1683b59SHong Zhang 
1000b28d80bdSHong Zhang   cm = c->m;
1001b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1002*e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1003b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1004b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1005b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1006b1683b59SHong Zhang     c->ncolumns[i] = n;
1007b1683b59SHong Zhang     if (n) {
1008d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1009b1683b59SHong Zhang     }
1010d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1011d2853540SHong Zhang     columns_i       += n;
1012b1683b59SHong Zhang 
1013b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1014e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1015*e99be685SHong Zhang 
1016b1683b59SHong Zhang     /* loop over columns*/
1017b1683b59SHong Zhang     for (j=0; j<n; j++) {
1018b1683b59SHong Zhang       col     = is[j];
1019d2853540SHong Zhang       row_idx = cj + ci[col];
1020b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1021b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1022b1683b59SHong Zhang       for (k=0; k<m; k++) {
1023*e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1024d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1025b1683b59SHong Zhang       }
1026b1683b59SHong Zhang     }
1027b1683b59SHong Zhang     /* count the number of hits */
1028b1683b59SHong Zhang     nrows = 0;
1029e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1030b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1031b1683b59SHong Zhang     }
1032b1683b59SHong Zhang     c->nrows[i]      = nrows;
1033d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1034d2853540SHong Zhang 
1035b1683b59SHong Zhang     nrows       = 0;
1036e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1037b1683b59SHong Zhang       if (rowhit[j]) {
1038d2853540SHong Zhang         rows_i[nrows]            = j;
1039*e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1040b1683b59SHong Zhang         nrows++;
1041b1683b59SHong Zhang       }
1042b1683b59SHong Zhang     }
1043b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1044d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1045b1683b59SHong Zhang   }
1046*e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1047b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1048b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1049d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1050d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1051d2853540SHong Zhang #endif
1052b28d80bdSHong Zhang 
1053d2853540SHong Zhang   c->colorforrow     = colorforrow;
1054d2853540SHong Zhang   c->rows            = rows;
1055d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1056d2853540SHong Zhang   c->colorforcol     = colorforcol;
1057d2853540SHong Zhang   c->columns         = columns;
1058f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1059b1683b59SHong Zhang   PetscFunctionReturn(0);
1060b1683b59SHong Zhang }
1061