xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
282412ba7SHong Zhang 
342c57489SHong Zhang /*
442c57489SHong Zhang   Defines projective product routines where A is a MPIAIJ matrix
542c57489SHong Zhang           C = P^T * A * P
642c57489SHong Zhang */
742c57489SHong Zhang 
87c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
97c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
107c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
1142c57489SHong Zhang #include "petscbt.h"
1242c57489SHong Zhang 
13*09573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
14cc6cb787SHong Zhang #undef __FUNCT__
15cc6cb787SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatPtAP"
16cc6cb787SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_MatPtAP(Mat A)
17cc6cb787SHong Zhang {
18cc6cb787SHong Zhang   PetscErrorCode       ierr;
19cc6cb787SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
20776b82aeSLisandro Dalcin   PetscContainer       container;
21cc6cb787SHong Zhang 
22cc6cb787SHong Zhang   PetscFunctionBegin;
23cc6cb787SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
24cc6cb787SHong Zhang   if (container) {
25776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
26cc6cb787SHong Zhang     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
27cc6cb787SHong Zhang     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
28cc6cb787SHong Zhang     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
29cc6cb787SHong Zhang     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
30cc6cb787SHong Zhang     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
31c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
32cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
33c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
34cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
3505b42c5fSBarry Smith     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
3605b42c5fSBarry Smith     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
3705b42c5fSBarry Smith     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
3826283091SBarry Smith     ierr = PetscLayoutDestroy(merge->rowmap);CHKERRQ(ierr);
39cc6cb787SHong Zhang 
40776b82aeSLisandro Dalcin     ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
41cc6cb787SHong Zhang     ierr = PetscObjectCompose((PetscObject)A,"MatMergeSeqsToMPI",0);CHKERRQ(ierr);
42cc6cb787SHong Zhang   }
43cc6cb787SHong Zhang   ierr = merge->MatDestroy(A);CHKERRQ(ierr);
44cc6cb787SHong Zhang   ierr = PetscFree(merge);CHKERRQ(ierr);
45cc6cb787SHong Zhang   PetscFunctionReturn(0);
46cc6cb787SHong Zhang }
47cc6cb787SHong Zhang 
48cc6cb787SHong Zhang #undef __FUNCT__
49cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
50f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
51f0c0a3a6SBarry Smith {
52cc6cb787SHong Zhang   PetscErrorCode       ierr;
53cc6cb787SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
54776b82aeSLisandro Dalcin   PetscContainer       container;
55cc6cb787SHong Zhang 
56cc6cb787SHong Zhang   PetscFunctionBegin;
57cc6cb787SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
58cc6cb787SHong Zhang   if (container) {
59776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
60e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
61cc6cb787SHong Zhang   ierr = (*merge->MatDuplicate)(A,op,M);CHKERRQ(ierr);
62cc6cb787SHong Zhang   (*M)->ops->destroy   = merge->MatDestroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
63cc6cb787SHong Zhang   (*M)->ops->duplicate = merge->MatDuplicate; /* =MatDuplicate_ MPIAIJ */
64cc6cb787SHong Zhang   PetscFunctionReturn(0);
65cc6cb787SHong Zhang }
66cc6cb787SHong Zhang 
6742c57489SHong Zhang #undef __FUNCT__
687a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
697a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
7042c57489SHong Zhang {
7142c57489SHong Zhang   PetscErrorCode ierr;
7242c57489SHong Zhang 
7342c57489SHong Zhang   PetscFunctionBegin;
7465e19b50SBarry Smith   if (!P->ops->ptapsymbolic_mpiaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
757a7894deSKris Buschelman   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
767a7894deSKris Buschelman   PetscFunctionReturn(0);
777a7894deSKris Buschelman }
787a7894deSKris Buschelman 
797a7894deSKris Buschelman #undef __FUNCT__
807a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
817a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
827a7894deSKris Buschelman {
837a7894deSKris Buschelman   PetscErrorCode ierr;
847a7894deSKris Buschelman 
857a7894deSKris Buschelman   PetscFunctionBegin;
8665e19b50SBarry Smith   if (!P->ops->ptapnumeric_mpiaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
877a7894deSKris Buschelman   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
8842c57489SHong Zhang   PetscFunctionReturn(0);
8942c57489SHong Zhang }
9042c57489SHong Zhang 
9142c57489SHong Zhang #undef __FUNCT__
9242c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9342c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9442c57489SHong Zhang {
9542c57489SHong Zhang   PetscErrorCode       ierr;
966c6e00e0SHong Zhang   Mat                  B_mpi;
97de0260b3SHong Zhang   Mat_MatMatMultMPI    *ap;
98776b82aeSLisandro Dalcin   PetscContainer       container;
99a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
10083445d95SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
10142c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10242c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
103ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10483445d95SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
10513f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
106d0f46423SBarry Smith   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
10742c57489SHong Zhang   PetscBT              lnkbt;
1087adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
109ded4f561SHong Zhang   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
11042c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11142c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
112ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11342c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
114ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
115ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11642c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
117f2b054eeSHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0;
11813f74950SBarry Smith   PetscMPIInt          j;
11942c57489SHong Zhang 
12042c57489SHong Zhang   PetscFunctionBegin;
12183445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
12283445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12383445d95SHong Zhang 
1246c6e00e0SHong Zhang   /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */
1256c6e00e0SHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
1266c6e00e0SHong Zhang   if (container) {
1274067f9b5SMatthew Knepley     /* reset functions */
1284067f9b5SMatthew Knepley     ierr = PetscContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr);
1294067f9b5SMatthew Knepley     P->ops->destroy = ap->MatDestroy;
1304067f9b5SMatthew Knepley     P->ops->duplicate = ap->MatDuplicate;
1314067f9b5SMatthew Knepley     /* destroy container and contents */
132776b82aeSLisandro Dalcin     ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
1336c6e00e0SHong Zhang     ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
1346c6e00e0SHong Zhang   }
1356c6e00e0SHong Zhang 
1366c6e00e0SHong Zhang   /* create the container 'Mat_MatMatMultMPI' and attach it to P */
1376c6e00e0SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr);
1386c6e00e0SHong Zhang   ap->abi=PETSC_NULL; ap->abj=PETSC_NULL;
1396c6e00e0SHong Zhang   ap->abnz_max = 0;
1406c6e00e0SHong Zhang 
141776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
142776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,ap);CHKERRQ(ierr);
1436c6e00e0SHong Zhang   ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
1446c6e00e0SHong Zhang   ap->MatDestroy  = P->ops->destroy;
1456c6e00e0SHong Zhang   P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
1466c6e00e0SHong Zhang   ap->reuse       = MAT_INITIAL_MATRIX;
147776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
1486c6e00e0SHong Zhang 
1496c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1501d79065fSBarry Smith   ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
1516c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
1526c6e00e0SHong Zhang   ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr);
1536c6e00e0SHong Zhang 
1546c6e00e0SHong Zhang   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
1556c6e00e0SHong Zhang   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
1566c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1576c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
15842c57489SHong Zhang 
159fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
160fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
16182412ba7SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
16282412ba7SHong Zhang   ap->abi = api;
16382412ba7SHong Zhang   api[0] = 0;
16483445d95SHong Zhang 
16583445d95SHong Zhang   /* create and initialize a linked list */
16683445d95SHong Zhang   nlnk = pN+1;
16783445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
16883445d95SHong Zhang 
16982412ba7SHong Zhang   /* Initial FreeSpace size is fill*nnz(A) */
170a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
171f4a743e1SHong Zhang   current_space = free_space;
172f4a743e1SHong Zhang 
173f4a743e1SHong Zhang   for (i=0;i<am;i++) {
17482412ba7SHong Zhang     apnz = 0;
175f4a743e1SHong Zhang     /* diagonal portion of A */
176ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
177ded4f561SHong Zhang     for (j=0; j<nzi; j++){
17882412ba7SHong Zhang       row = *adj++;
17982412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
180ded4f561SHong Zhang       Jptr  = pj_loc + pi_loc[row];
18183445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
182ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
18382412ba7SHong Zhang       apnz += nlnk;
184f4a743e1SHong Zhang     }
185f4a743e1SHong Zhang     /* off-diagonal portion of A */
186ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
187ded4f561SHong Zhang     for (j=0; j<nzi; j++){
18882412ba7SHong Zhang       row = *aoj++;
18982412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
190ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
191ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
19282412ba7SHong Zhang       apnz += nlnk;
193f4a743e1SHong Zhang     }
194f4a743e1SHong Zhang 
19582412ba7SHong Zhang     api[i+1] = api[i] + apnz;
19682412ba7SHong Zhang     if (ap->abnz_max < apnz) ap->abnz_max = apnz;
197f4a743e1SHong Zhang 
19883445d95SHong Zhang     /* if free space is not available, double the total space in the list */
19982412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
2002ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
201f2b054eeSHong Zhang       nspacedouble++;
202f4a743e1SHong Zhang     }
203f4a743e1SHong Zhang 
204f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
20582412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
20682412ba7SHong Zhang     current_space->array           += apnz;
20782412ba7SHong Zhang     current_space->local_used      += apnz;
20882412ba7SHong Zhang     current_space->local_remaining -= apnz;
209f4a743e1SHong Zhang   }
21082412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
211f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
21282412ba7SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr);
21382412ba7SHong Zhang   apj = ap->abj;
214a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr);
215f4a743e1SHong Zhang 
216ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
217ded4f561SHong Zhang   /*----------------------------------------------------*/
218de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
21942c57489SHong Zhang 
220ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
221d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
22283445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
223de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
224de0260b3SHong Zhang   coi[0] = 0;
22542c57489SHong Zhang 
22682412ba7SHong Zhang   /* set initial free space to be 3*pon*max( nnz(AP) per row) */
22782412ba7SHong Zhang   nnz           = 3*pon*ap->abnz_max + 1;
228a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
22942c57489SHong Zhang   current_space = free_space;
23042c57489SHong Zhang 
231de0260b3SHong Zhang   for (i=0; i<pon; i++) {
232ded4f561SHong Zhang     nnz  = 0;
23382412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
23482412ba7SHong Zhang     j     = pnz;
235de0260b3SHong Zhang     ptJ   = potj + poti[i+1];
23683445d95SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
23783445d95SHong Zhang       j--; ptJ--;
23882412ba7SHong Zhang       row  = *ptJ; /* row of AP == col of Pot */
23982412ba7SHong Zhang       apnz = api[row+1] - api[row];
240ded4f561SHong Zhang       Jptr   = apj + api[row];
24183445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
242ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
243ded4f561SHong Zhang       nnz += nlnk;
24442c57489SHong Zhang     }
24542c57489SHong Zhang 
24683445d95SHong Zhang     /* If free space is not available, double the total space in the list */
247ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2484238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
24942c57489SHong Zhang     }
25042c57489SHong Zhang 
25142c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
252ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
253ded4f561SHong Zhang     current_space->array           += nnz;
254ded4f561SHong Zhang     current_space->local_used      += nnz;
255ded4f561SHong Zhang     current_space->local_remaining -= nnz;
256ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
25742c57489SHong Zhang   }
258de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
259a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
260de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
26142c57489SHong Zhang 
262ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
263ded4f561SHong Zhang   /*----------------------------------------------*/
26442c57489SHong Zhang   /* determine row ownership */
26583445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
26626283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2677a2fc3feSBarry Smith   merge->rowmap->n = pn;
2687a2fc3feSBarry Smith   merge->rowmap->bs = 1;
26926283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2707a2fc3feSBarry Smith   owners = merge->rowmap->range;
27142c57489SHong Zhang 
27242c57489SHong Zhang   /* determine the number of messages to send, their lengths */
27342c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
27483445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
27542c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
27642c57489SHong Zhang   len_s = merge->len_s;
27742c57489SHong Zhang   merge->nsend = 0;
27883445d95SHong Zhang 
279de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
280de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
281de0260b3SHong Zhang 
28283445d95SHong Zhang   proc = 0;
283de0260b3SHong Zhang   for (i=0; i<pon; i++){
284de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
285de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
286de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
28783445d95SHong Zhang   }
288de0260b3SHong Zhang 
289de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
290de0260b3SHong Zhang   owners_co[0] = 0;
29142c57489SHong Zhang   for (proc=0; proc<size; proc++){
292de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
29383445d95SHong Zhang     if (len_si[proc]){
29442c57489SHong Zhang       merge->nsend++;
29583445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
29642c57489SHong Zhang       len += len_si[proc];
29742c57489SHong Zhang     }
29842c57489SHong Zhang   }
29942c57489SHong Zhang 
300ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
30142c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
30242c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
30342c57489SHong Zhang 
304ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
305fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
306ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
30742c57489SHong Zhang 
308ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
30942c57489SHong Zhang 
31042c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
31142c57489SHong Zhang     if (!len_s[proc]) continue;
312de0260b3SHong Zhang     i = owners_co[proc];
313ded4f561SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
31442c57489SHong Zhang     k++;
31542c57489SHong Zhang   }
31642c57489SHong Zhang 
317ded4f561SHong Zhang   /* receives and sends of coj are complete */
318ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
319ded4f561SHong Zhang   i = merge->nrecv;
320ded4f561SHong Zhang   while (i--) {
321ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
322ded4f561SHong Zhang   }
323ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3240c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
32582412ba7SHong Zhang 
326ded4f561SHong Zhang   /* send and recv coi */
327ded4f561SHong Zhang   /*-------------------*/
328ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
32942c57489SHong Zhang 
33042c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
33142c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
33242c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
33342c57489SHong Zhang     if (!len_s[proc]) continue;
33442c57489SHong Zhang     /* form outgoing message for i-structure:
33542c57489SHong Zhang          buf_si[0]:                 nrows to be sent
33642c57489SHong Zhang                [1:nrows]:           row index (global)
33742c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
33842c57489SHong Zhang     */
33942c57489SHong Zhang     /*-------------------------------------------*/
34042c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
34142c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
34242c57489SHong Zhang     buf_si[0]   = nrows;
34342c57489SHong Zhang     buf_si_i[0] = 0;
34442c57489SHong Zhang     nrows = 0;
345de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
346ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
347ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
348de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
34942c57489SHong Zhang       nrows++;
35042c57489SHong Zhang     }
351ded4f561SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
35242c57489SHong Zhang     k++;
35342c57489SHong Zhang     buf_si += len_si[proc];
35442c57489SHong Zhang   }
355ded4f561SHong Zhang   i = merge->nrecv;
356ded4f561SHong Zhang   while (i--) {
357ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
358ded4f561SHong Zhang   }
359ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3600c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
361f2b054eeSHong Zhang   /*
362ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
36342c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
364ae15b995SBarry Smith     ierr = PetscInfo3(A,"recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr);
36542c57489SHong Zhang   }
366f2b054eeSHong Zhang   */
36742c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
36842c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
369ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
370ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
37142c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
37242c57489SHong Zhang 
373de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
374de0260b3SHong Zhang   /*------------------------------------------*/
375ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
376ded4f561SHong Zhang 
37742c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
37842c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
37942c57489SHong Zhang   bi[0] = 0;
38042c57489SHong Zhang 
381ded4f561SHong Zhang   /* set initial free space to be 3*pn*max( nnz(AP) per row) */
382ded4f561SHong Zhang   nnz           = 3*pn*ap->abnz_max + 1;
383a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
38442c57489SHong Zhang   current_space = free_space;
38542c57489SHong Zhang 
386687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
38742c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
38842c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
38942c57489SHong Zhang     nrows       = *buf_ri_k[k];
39042c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
39142c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
39242c57489SHong Zhang   }
39342c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
39442c57489SHong Zhang   for (i=0; i<pn; i++) {
395ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
396ded4f561SHong Zhang     nnz = 0;
397ded4f561SHong Zhang     pnz  = pdti[i+1] - pdti[i];
398ded4f561SHong Zhang     j    = pnz;
399ded4f561SHong Zhang     ptJ  = pdtj + pdti[i+1];
400ded4f561SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
401ded4f561SHong Zhang       j--; ptJ--;
402ded4f561SHong Zhang       row  = *ptJ; /* row of AP == col of Pt */
403ded4f561SHong Zhang       apnz = api[row+1] - api[row];
404ded4f561SHong Zhang       Jptr   = apj + api[row];
405ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
406ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
407ded4f561SHong Zhang       nnz += nlnk;
408ded4f561SHong Zhang     }
40942c57489SHong Zhang     /* add received col data into lnk */
41042c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
41142c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
412ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
413ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
414ded4f561SHong Zhang         ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
415ded4f561SHong Zhang         nnz += nlnk;
41642c57489SHong Zhang         nextrow[k]++; nextci[k]++;
41742c57489SHong Zhang       }
41842c57489SHong Zhang     }
41942c57489SHong Zhang 
42042c57489SHong Zhang     /* if free space is not available, make more free space */
421ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4224238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
42342c57489SHong Zhang     }
42442c57489SHong Zhang     /* copy data into free space, then initialize lnk */
425ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
426ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
427ded4f561SHong Zhang     current_space->array           += nnz;
428ded4f561SHong Zhang     current_space->local_used      += nnz;
429ded4f561SHong Zhang     current_space->local_remaining -= nnz;
430ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
43142c57489SHong Zhang   }
432ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4330572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
43442c57489SHong Zhang 
43542c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
436a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
43742c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
43842c57489SHong Zhang 
43942c57489SHong Zhang   /* create symbolic parallel matrix B_mpi */
44042c57489SHong Zhang   /*---------------------------------------*/
441f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
442f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
44342c57489SHong Zhang   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
44442c57489SHong Zhang   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
44542c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
44642c57489SHong Zhang 
44742c57489SHong Zhang   merge->bi            = bi;
44842c57489SHong Zhang   merge->bj            = bj;
449de0260b3SHong Zhang   merge->coi           = coi;
450de0260b3SHong Zhang   merge->coj           = coj;
45142c57489SHong Zhang   merge->buf_ri        = buf_ri;
45242c57489SHong Zhang   merge->buf_rj        = buf_rj;
453de0260b3SHong Zhang   merge->owners_co     = owners_co;
454cc6cb787SHong Zhang   merge->MatDestroy    = B_mpi->ops->destroy;
455cc6cb787SHong Zhang   merge->MatDuplicate  = B_mpi->ops->duplicate;
456cc6cb787SHong Zhang 
457cc6cb787SHong Zhang   /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
458cc6cb787SHong Zhang   B_mpi->assembled     = PETSC_FALSE;
459cc6cb787SHong Zhang   B_mpi->ops->destroy  = MatDestroy_MPIAIJ_MatPtAP;
460cc6cb787SHong Zhang   B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
46142c57489SHong Zhang 
46242c57489SHong Zhang   /* attach the supporting struct to B_mpi for reuse */
463776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
464776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,merge);CHKERRQ(ierr);
46542c57489SHong Zhang   ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr);
46642c57489SHong Zhang   *C = B_mpi;
467f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
468f2b054eeSHong Zhang   if (bi[pn] != 0) {
469f2b054eeSHong Zhang     PetscReal afill = ((PetscReal)bi[pn])/(adi[am]+aoi[am]);
470f2b054eeSHong Zhang     if (afill < 1.0) afill = 1.0;
471f2b054eeSHong Zhang     ierr = PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G when computing A*P.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
472f2b054eeSHong Zhang     ierr = PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
473f2b054eeSHong Zhang   } else {
474f2b054eeSHong Zhang     ierr = PetscInfo(B_mpi,"Empty matrix product\n");CHKERRQ(ierr);
475f2b054eeSHong Zhang   }
476f2b054eeSHong Zhang #endif
47742c57489SHong Zhang   PetscFunctionReturn(0);
47842c57489SHong Zhang }
47942c57489SHong Zhang 
48042c57489SHong Zhang #undef __FUNCT__
48142c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
48242c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
48342c57489SHong Zhang {
48442c57489SHong Zhang   PetscErrorCode       ierr;
48542c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
486de0260b3SHong Zhang   Mat_MatMatMultMPI    *ap;
487776b82aeSLisandro Dalcin   PetscContainer       cont_merge,cont_ptap;
488de0260b3SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
48942c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
490de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
49142c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
492bf3909cdSBarry Smith   PetscInt             *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp;
49382412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
49482412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
49582412ba7SHong Zhang   MatScalar            *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth;
496d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
4977adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
49842c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
49982412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
50042c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
50150f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
50242c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
50342c57489SHong Zhang   MPI_Status           *status;
50482412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
50582412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
506d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
50742c57489SHong Zhang 
50842c57489SHong Zhang   PetscFunctionBegin;
50942c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
51042c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
51142c57489SHong Zhang 
51242c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr);
51342c57489SHong Zhang   if (cont_merge) {
514776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr);
515e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
5166c6e00e0SHong Zhang 
517f33d1a9aSHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr);
51842c57489SHong Zhang   if (cont_ptap) {
519776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr);
5206c6e00e0SHong Zhang     if (ap->reuse == MAT_INITIAL_MATRIX){
5216c6e00e0SHong Zhang       ap->reuse = MAT_REUSE_MATRIX;
5226c6e00e0SHong Zhang     } else { /* update numerical values of P_oth and P_loc */
5231d79065fSBarry Smith       ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
5246c6e00e0SHong Zhang       ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr);
5256c6e00e0SHong Zhang     }
526e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
52742c57489SHong Zhang 
52842c57489SHong Zhang   /* get data from symbolic products */
529de0260b3SHong Zhang   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
530de0260b3SHong Zhang   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
53182412ba7SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc;
53242c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
53342c57489SHong Zhang 
534de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
535de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
536de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
537de0260b3SHong Zhang 
538de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5397a2fc3feSBarry Smith   owners = merge->rowmap->range;
540de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
541de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
54242c57489SHong Zhang 
543f4a743e1SHong Zhang   /* get data from symbolic A*P */
544de0260b3SHong Zhang   ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
545de0260b3SHong Zhang   ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
54642c57489SHong Zhang 
547f4a743e1SHong Zhang   /* compute numeric C_seq=P_loc^T*A_loc*P */
54882412ba7SHong Zhang   api = ap->abi; apj = ap->abj;
54942c57489SHong Zhang   for (i=0;i<am;i++) {
550f4a743e1SHong Zhang     /* form i-th sparse row of A*P */
55182412ba7SHong Zhang     apnz = api[i+1] - api[i];
55282412ba7SHong Zhang     apJ  = apj + api[i];
55342c57489SHong Zhang     /* diagonal portion of A */
55482412ba7SHong Zhang     anz  = adi[i+1] - adi[i];
55582412ba7SHong Zhang     for (j=0;j<anz;j++) {
55682412ba7SHong Zhang       row = *adj++;
55782412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
55882412ba7SHong Zhang       pj  = pj_loc + pi_loc[row];
55982412ba7SHong Zhang       pa  = pa_loc + pi_loc[row];
560f4a743e1SHong Zhang       nextp = 0;
56182412ba7SHong Zhang       for (k=0; nextp<pnz; k++) {
56282412ba7SHong Zhang         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
56382412ba7SHong Zhang           apa[k] += (*ada)*pa[nextp++];
56442c57489SHong Zhang         }
56542c57489SHong Zhang       }
566dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
56742c57489SHong Zhang       ada++;
56842c57489SHong Zhang     }
56942c57489SHong Zhang     /* off-diagonal portion of A */
57082412ba7SHong Zhang     anz  = aoi[i+1] - aoi[i];
57182412ba7SHong Zhang     for (j=0; j<anz; j++) {
57282412ba7SHong Zhang       row = *aoj++;
57382412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
57482412ba7SHong Zhang       pj  = pj_oth + pi_oth[row];
57582412ba7SHong Zhang       pa  = pa_oth + pi_oth[row];
576f4a743e1SHong Zhang       nextp = 0;
57782412ba7SHong Zhang       for (k=0; nextp<pnz; k++) {
57882412ba7SHong Zhang         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
57982412ba7SHong Zhang           apa[k] += (*aoa)*pa[nextp++];
58042c57489SHong Zhang         }
58142c57489SHong Zhang       }
582dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
58342c57489SHong Zhang       aoa++;
58442c57489SHong Zhang     }
58542c57489SHong Zhang 
58682412ba7SHong Zhang     /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
58782412ba7SHong Zhang     pnz = pi_loc[i+1] - pi_loc[i];
58882412ba7SHong Zhang     for (j=0; j<pnz; j++) {
58942c57489SHong Zhang       nextap = 0;
59082412ba7SHong Zhang       row    = *pJ++; /* global index */
59182412ba7SHong Zhang       if (row < pcstart || row >=pcend) { /* put the value into Co */
59282412ba7SHong Zhang         cj  = coj + coi[*poJ];
59382412ba7SHong Zhang         ca  = coa + coi[*poJ++];
59482412ba7SHong Zhang       } else {                            /* put the value into Cd */
59582412ba7SHong Zhang         cj   = bj + bi[*pdJ];
59682412ba7SHong Zhang         ca   = ba + bi[*pdJ++];
59742c57489SHong Zhang       }
59882412ba7SHong Zhang       for (k=0; nextap<apnz; k++) {
59982412ba7SHong Zhang         if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++];
60042c57489SHong Zhang       }
601dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
60282412ba7SHong Zhang       pA++;
603de0260b3SHong Zhang     }
604de0260b3SHong Zhang 
60542c57489SHong Zhang     /* zero the current row info for A*P */
60682412ba7SHong Zhang     ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
60742c57489SHong Zhang   }
60842c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
60942c57489SHong Zhang 
610de0260b3SHong Zhang   /* send and recv matrix values */
611de0260b3SHong Zhang   /*-----------------------------*/
61242c57489SHong Zhang   buf_ri = merge->buf_ri;
61342c57489SHong Zhang   buf_rj = merge->buf_rj;
61442c57489SHong Zhang   len_s  = merge->len_s;
615fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
61642c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
61742c57489SHong Zhang 
61842c57489SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr);
61942c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
62042c57489SHong Zhang     if (!len_s[proc]) continue;
621de0260b3SHong Zhang     i = merge->owners_co[proc];
622de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
62342c57489SHong Zhang     k++;
62442c57489SHong Zhang   }
62582412ba7SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
6260c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
6270c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
62842c57489SHong Zhang   ierr = PetscFree(status);CHKERRQ(ierr);
62942c57489SHong Zhang 
63042c57489SHong Zhang   ierr = PetscFree(s_waits);CHKERRQ(ierr);
63142c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
63282412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
63342c57489SHong Zhang 
63482412ba7SHong Zhang   /* insert local and received values into C */
63582412ba7SHong Zhang   /*-----------------------------------------*/
636687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
63742c57489SHong Zhang 
63842c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
63942c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
64042c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
64142c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
64282412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
64342c57489SHong Zhang   }
64442c57489SHong Zhang 
645de0260b3SHong Zhang   for (i=0; i<cm; i++) {
64682412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
64742c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
648de0260b3SHong Zhang     ba_i = ba + bi[i];
64982412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
65042c57489SHong Zhang     /* add received vals into ba */
65142c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
65242c57489SHong Zhang       /* i-th row */
65342c57489SHong Zhang       if (i == *nextrow[k]) {
65482412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
65582412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
65682412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
65782412ba7SHong Zhang         nextcj = 0;
65882412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
65982412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
66082412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
66142c57489SHong Zhang           }
66242c57489SHong Zhang         }
66382412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
66442c57489SHong Zhang       }
66542c57489SHong Zhang     }
66682412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
667dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
66842c57489SHong Zhang   }
669fd0ff01cSHong Zhang   ierr = MatSetBlockSize(C,1);CHKERRQ(ierr);
67042c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67142c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67242c57489SHong Zhang 
673de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
674c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
67542c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
6760572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
67742c57489SHong Zhang   PetscFunctionReturn(0);
67842c57489SHong Zhang }
679