xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision c05d87d6c2fa1d0952fe676de20339d09daa79bc)
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 
13cc6cb787SHong Zhang 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);
31*c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
32cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
33*c05d87d6SBarry 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);
60cc6cb787SHong Zhang   } else {
61cc6cb787SHong Zhang     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
62cc6cb787SHong Zhang   }
63cc6cb787SHong Zhang   ierr = (*merge->MatDuplicate)(A,op,M);CHKERRQ(ierr);
64cc6cb787SHong Zhang   (*M)->ops->destroy   = merge->MatDestroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
65cc6cb787SHong Zhang   (*M)->ops->duplicate = merge->MatDuplicate; /* =MatDuplicate_ MPIAIJ */
66cc6cb787SHong Zhang   PetscFunctionReturn(0);
67cc6cb787SHong Zhang }
68cc6cb787SHong Zhang 
6942c57489SHong Zhang #undef __FUNCT__
707a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
717a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
7242c57489SHong Zhang {
7342c57489SHong Zhang   PetscErrorCode ierr;
7442c57489SHong Zhang 
7542c57489SHong Zhang   PetscFunctionBegin;
767a7894deSKris Buschelman   if (!P->ops->ptapsymbolic_mpiaij) {
777adad957SLisandro Dalcin     SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
7842c57489SHong Zhang   }
797a7894deSKris Buschelman   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
807a7894deSKris Buschelman   PetscFunctionReturn(0);
817a7894deSKris Buschelman }
827a7894deSKris Buschelman 
837a7894deSKris Buschelman #undef __FUNCT__
847a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
857a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
867a7894deSKris Buschelman {
877a7894deSKris Buschelman   PetscErrorCode ierr;
887a7894deSKris Buschelman 
897a7894deSKris Buschelman   PetscFunctionBegin;
907a7894deSKris Buschelman   if (!P->ops->ptapnumeric_mpiaij) {
917adad957SLisandro Dalcin     SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
927a7894deSKris Buschelman   }
937a7894deSKris Buschelman   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
9442c57489SHong Zhang   PetscFunctionReturn(0);
9542c57489SHong Zhang }
9642c57489SHong Zhang 
9742c57489SHong Zhang #undef __FUNCT__
9842c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9942c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
10042c57489SHong Zhang {
10142c57489SHong Zhang   PetscErrorCode       ierr;
1026c6e00e0SHong Zhang   Mat                  B_mpi;
103de0260b3SHong Zhang   Mat_MatMatMultMPI    *ap;
104776b82aeSLisandro Dalcin   PetscContainer       container;
105a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
10683445d95SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
10742c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10842c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
109ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
11083445d95SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
11113f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
112d0f46423SBarry Smith   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
11342c57489SHong Zhang   PetscBT              lnkbt;
1147adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
115ded4f561SHong Zhang   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
11642c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11742c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
118ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11942c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
120ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
121ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
12242c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
123f2b054eeSHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0;
12413f74950SBarry Smith   PetscMPIInt          j;
12542c57489SHong Zhang 
12642c57489SHong Zhang   PetscFunctionBegin;
12783445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
12883445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12983445d95SHong Zhang 
1306c6e00e0SHong Zhang   /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */
1316c6e00e0SHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
1326c6e00e0SHong Zhang   if (container) {
1334067f9b5SMatthew Knepley     /* reset functions */
1344067f9b5SMatthew Knepley     ierr = PetscContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr);
1354067f9b5SMatthew Knepley     P->ops->destroy = ap->MatDestroy;
1364067f9b5SMatthew Knepley     P->ops->duplicate = ap->MatDuplicate;
1374067f9b5SMatthew Knepley     /* destroy container and contents */
138776b82aeSLisandro Dalcin     ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
1396c6e00e0SHong Zhang     ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
1406c6e00e0SHong Zhang   }
1416c6e00e0SHong Zhang 
1426c6e00e0SHong Zhang   /* create the container 'Mat_MatMatMultMPI' and attach it to P */
1436c6e00e0SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr);
1446c6e00e0SHong Zhang   ap->abi=PETSC_NULL; ap->abj=PETSC_NULL;
1456c6e00e0SHong Zhang   ap->abnz_max = 0;
1466c6e00e0SHong Zhang 
147776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
148776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,ap);CHKERRQ(ierr);
1496c6e00e0SHong Zhang   ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
1506c6e00e0SHong Zhang   ap->MatDestroy  = P->ops->destroy;
1516c6e00e0SHong Zhang   P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
1526c6e00e0SHong Zhang   ap->reuse       = MAT_INITIAL_MATRIX;
153776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
1546c6e00e0SHong Zhang 
1556c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1561d79065fSBarry Smith   ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
1576c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
1586c6e00e0SHong Zhang   ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr);
1596c6e00e0SHong Zhang 
1606c6e00e0SHong Zhang   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
1616c6e00e0SHong Zhang   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
1626c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1636c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
16442c57489SHong Zhang 
165fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
166fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
16782412ba7SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
16882412ba7SHong Zhang   ap->abi = api;
16982412ba7SHong Zhang   api[0] = 0;
17083445d95SHong Zhang 
17183445d95SHong Zhang   /* create and initialize a linked list */
17283445d95SHong Zhang   nlnk = pN+1;
17383445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
17483445d95SHong Zhang 
17582412ba7SHong Zhang   /* Initial FreeSpace size is fill*nnz(A) */
176a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
177f4a743e1SHong Zhang   current_space = free_space;
178f4a743e1SHong Zhang 
179f4a743e1SHong Zhang   for (i=0;i<am;i++) {
18082412ba7SHong Zhang     apnz = 0;
181f4a743e1SHong Zhang     /* diagonal portion of A */
182ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
183ded4f561SHong Zhang     for (j=0; j<nzi; j++){
18482412ba7SHong Zhang       row = *adj++;
18582412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
186ded4f561SHong Zhang       Jptr  = pj_loc + pi_loc[row];
18783445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
188ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
18982412ba7SHong Zhang       apnz += nlnk;
190f4a743e1SHong Zhang     }
191f4a743e1SHong Zhang     /* off-diagonal portion of A */
192ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
193ded4f561SHong Zhang     for (j=0; j<nzi; j++){
19482412ba7SHong Zhang       row = *aoj++;
19582412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
196ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
197ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
19882412ba7SHong Zhang       apnz += nlnk;
199f4a743e1SHong Zhang     }
200f4a743e1SHong Zhang 
20182412ba7SHong Zhang     api[i+1] = api[i] + apnz;
20282412ba7SHong Zhang     if (ap->abnz_max < apnz) ap->abnz_max = apnz;
203f4a743e1SHong Zhang 
20483445d95SHong Zhang     /* if free space is not available, double the total space in the list */
20582412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
2062ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
207f2b054eeSHong Zhang       nspacedouble++;
208f4a743e1SHong Zhang     }
209f4a743e1SHong Zhang 
210f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
21182412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
21282412ba7SHong Zhang     current_space->array           += apnz;
21382412ba7SHong Zhang     current_space->local_used      += apnz;
21482412ba7SHong Zhang     current_space->local_remaining -= apnz;
215f4a743e1SHong Zhang   }
21682412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
217f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
21882412ba7SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr);
21982412ba7SHong Zhang   apj = ap->abj;
220a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr);
221f4a743e1SHong Zhang 
222ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
223ded4f561SHong Zhang   /*----------------------------------------------------*/
224de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
22542c57489SHong Zhang 
226ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
227d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
22883445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
229de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
230de0260b3SHong Zhang   coi[0] = 0;
23142c57489SHong Zhang 
23282412ba7SHong Zhang   /* set initial free space to be 3*pon*max( nnz(AP) per row) */
23382412ba7SHong Zhang   nnz           = 3*pon*ap->abnz_max + 1;
234a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
23542c57489SHong Zhang   current_space = free_space;
23642c57489SHong Zhang 
237de0260b3SHong Zhang   for (i=0; i<pon; i++) {
238ded4f561SHong Zhang     nnz  = 0;
23982412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
24082412ba7SHong Zhang     j     = pnz;
241de0260b3SHong Zhang     ptJ   = potj + poti[i+1];
24283445d95SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
24383445d95SHong Zhang       j--; ptJ--;
24482412ba7SHong Zhang       row  = *ptJ; /* row of AP == col of Pot */
24582412ba7SHong Zhang       apnz = api[row+1] - api[row];
246ded4f561SHong Zhang       Jptr   = apj + api[row];
24783445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
248ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
249ded4f561SHong Zhang       nnz += nlnk;
25042c57489SHong Zhang     }
25142c57489SHong Zhang 
25283445d95SHong Zhang     /* If free space is not available, double the total space in the list */
253ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2544238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
25542c57489SHong Zhang     }
25642c57489SHong Zhang 
25742c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
258ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
259ded4f561SHong Zhang     current_space->array           += nnz;
260ded4f561SHong Zhang     current_space->local_used      += nnz;
261ded4f561SHong Zhang     current_space->local_remaining -= nnz;
262ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
26342c57489SHong Zhang   }
264de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
265a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
266de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
26742c57489SHong Zhang 
268ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
269ded4f561SHong Zhang   /*----------------------------------------------*/
27042c57489SHong Zhang   /* determine row ownership */
27183445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
27226283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2737a2fc3feSBarry Smith   merge->rowmap->n = pn;
2747a2fc3feSBarry Smith   merge->rowmap->bs = 1;
27526283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2767a2fc3feSBarry Smith   owners = merge->rowmap->range;
27742c57489SHong Zhang 
27842c57489SHong Zhang   /* determine the number of messages to send, their lengths */
27942c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
28083445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
28142c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
28242c57489SHong Zhang   len_s = merge->len_s;
28342c57489SHong Zhang   merge->nsend = 0;
28483445d95SHong Zhang 
285de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
286de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
287de0260b3SHong Zhang 
28883445d95SHong Zhang   proc = 0;
289de0260b3SHong Zhang   for (i=0; i<pon; i++){
290de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
291de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
292de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
29383445d95SHong Zhang   }
294de0260b3SHong Zhang 
295de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
296de0260b3SHong Zhang   owners_co[0] = 0;
29742c57489SHong Zhang   for (proc=0; proc<size; proc++){
298de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
29983445d95SHong Zhang     if (len_si[proc]){
30042c57489SHong Zhang       merge->nsend++;
30183445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
30242c57489SHong Zhang       len += len_si[proc];
30342c57489SHong Zhang     }
30442c57489SHong Zhang   }
30542c57489SHong Zhang 
306ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
30742c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
30842c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
30942c57489SHong Zhang 
310ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
311fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
312ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
31342c57489SHong Zhang 
314ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
31542c57489SHong Zhang 
31642c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
31742c57489SHong Zhang     if (!len_s[proc]) continue;
318de0260b3SHong Zhang     i = owners_co[proc];
319ded4f561SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
32042c57489SHong Zhang     k++;
32142c57489SHong Zhang   }
32242c57489SHong Zhang 
323ded4f561SHong Zhang   /* receives and sends of coj are complete */
324ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
325ded4f561SHong Zhang   i = merge->nrecv;
326ded4f561SHong Zhang   while (i--) {
327ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
328ded4f561SHong Zhang   }
329ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3300c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
33182412ba7SHong Zhang 
332ded4f561SHong Zhang   /* send and recv coi */
333ded4f561SHong Zhang   /*-------------------*/
334ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
33542c57489SHong Zhang 
33642c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
33742c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
33842c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
33942c57489SHong Zhang     if (!len_s[proc]) continue;
34042c57489SHong Zhang     /* form outgoing message for i-structure:
34142c57489SHong Zhang          buf_si[0]:                 nrows to be sent
34242c57489SHong Zhang                [1:nrows]:           row index (global)
34342c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
34442c57489SHong Zhang     */
34542c57489SHong Zhang     /*-------------------------------------------*/
34642c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
34742c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
34842c57489SHong Zhang     buf_si[0]   = nrows;
34942c57489SHong Zhang     buf_si_i[0] = 0;
35042c57489SHong Zhang     nrows = 0;
351de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
352ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
353ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
354de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
35542c57489SHong Zhang       nrows++;
35642c57489SHong Zhang     }
357ded4f561SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
35842c57489SHong Zhang     k++;
35942c57489SHong Zhang     buf_si += len_si[proc];
36042c57489SHong Zhang   }
361ded4f561SHong Zhang   i = merge->nrecv;
362ded4f561SHong Zhang   while (i--) {
363ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
364ded4f561SHong Zhang   }
365ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3660c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
367f2b054eeSHong Zhang   /*
368ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
36942c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
370ae15b995SBarry 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);
37142c57489SHong Zhang   }
372f2b054eeSHong Zhang   */
37342c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
37442c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
375ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
376ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
37742c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
37842c57489SHong Zhang 
379de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
380de0260b3SHong Zhang   /*------------------------------------------*/
381ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
382ded4f561SHong Zhang 
38342c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
38442c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
38542c57489SHong Zhang   bi[0] = 0;
38642c57489SHong Zhang 
387ded4f561SHong Zhang   /* set initial free space to be 3*pn*max( nnz(AP) per row) */
388ded4f561SHong Zhang   nnz           = 3*pn*ap->abnz_max + 1;
389a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
39042c57489SHong Zhang   current_space = free_space;
39142c57489SHong Zhang 
392687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
39342c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
39442c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
39542c57489SHong Zhang     nrows       = *buf_ri_k[k];
39642c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
39742c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
39842c57489SHong Zhang   }
39942c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
40042c57489SHong Zhang   for (i=0; i<pn; i++) {
401ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
402ded4f561SHong Zhang     nnz = 0;
403ded4f561SHong Zhang     pnz  = pdti[i+1] - pdti[i];
404ded4f561SHong Zhang     j    = pnz;
405ded4f561SHong Zhang     ptJ  = pdtj + pdti[i+1];
406ded4f561SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
407ded4f561SHong Zhang       j--; ptJ--;
408ded4f561SHong Zhang       row  = *ptJ; /* row of AP == col of Pt */
409ded4f561SHong Zhang       apnz = api[row+1] - api[row];
410ded4f561SHong Zhang       Jptr   = apj + api[row];
411ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
412ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
413ded4f561SHong Zhang       nnz += nlnk;
414ded4f561SHong Zhang     }
41542c57489SHong Zhang     /* add received col data into lnk */
41642c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
41742c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
418ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
419ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
420ded4f561SHong Zhang         ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
421ded4f561SHong Zhang         nnz += nlnk;
42242c57489SHong Zhang         nextrow[k]++; nextci[k]++;
42342c57489SHong Zhang       }
42442c57489SHong Zhang     }
42542c57489SHong Zhang 
42642c57489SHong Zhang     /* if free space is not available, make more free space */
427ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4284238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
42942c57489SHong Zhang     }
43042c57489SHong Zhang     /* copy data into free space, then initialize lnk */
431ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
432ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
433ded4f561SHong Zhang     current_space->array           += nnz;
434ded4f561SHong Zhang     current_space->local_used      += nnz;
435ded4f561SHong Zhang     current_space->local_remaining -= nnz;
436ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
43742c57489SHong Zhang   }
438ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4390572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
44042c57489SHong Zhang 
44142c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
442a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
44342c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
44442c57489SHong Zhang 
44542c57489SHong Zhang   /* create symbolic parallel matrix B_mpi */
44642c57489SHong Zhang   /*---------------------------------------*/
447f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
448f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
44942c57489SHong Zhang   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
45042c57489SHong Zhang   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
45142c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
45242c57489SHong Zhang 
45342c57489SHong Zhang   merge->bi            = bi;
45442c57489SHong Zhang   merge->bj            = bj;
455de0260b3SHong Zhang   merge->coi           = coi;
456de0260b3SHong Zhang   merge->coj           = coj;
45742c57489SHong Zhang   merge->buf_ri        = buf_ri;
45842c57489SHong Zhang   merge->buf_rj        = buf_rj;
459de0260b3SHong Zhang   merge->owners_co     = owners_co;
460cc6cb787SHong Zhang   merge->MatDestroy    = B_mpi->ops->destroy;
461cc6cb787SHong Zhang   merge->MatDuplicate  = B_mpi->ops->duplicate;
462cc6cb787SHong Zhang 
463cc6cb787SHong Zhang   /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
464cc6cb787SHong Zhang   B_mpi->assembled     = PETSC_FALSE;
465cc6cb787SHong Zhang   B_mpi->ops->destroy  = MatDestroy_MPIAIJ_MatPtAP;
466cc6cb787SHong Zhang   B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
46742c57489SHong Zhang 
46842c57489SHong Zhang   /* attach the supporting struct to B_mpi for reuse */
469776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
470776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,merge);CHKERRQ(ierr);
47142c57489SHong Zhang   ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr);
47242c57489SHong Zhang   *C = B_mpi;
473f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
474f2b054eeSHong Zhang   if (bi[pn] != 0) {
475f2b054eeSHong Zhang     PetscReal afill = ((PetscReal)bi[pn])/(adi[am]+aoi[am]);
476f2b054eeSHong Zhang     if (afill < 1.0) afill = 1.0;
477f2b054eeSHong Zhang     ierr = PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G when computing A*P.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
478f2b054eeSHong Zhang     ierr = PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
479f2b054eeSHong Zhang   } else {
480f2b054eeSHong Zhang     ierr = PetscInfo(B_mpi,"Empty matrix product\n");CHKERRQ(ierr);
481f2b054eeSHong Zhang   }
482f2b054eeSHong Zhang #endif
48342c57489SHong Zhang   PetscFunctionReturn(0);
48442c57489SHong Zhang }
48542c57489SHong Zhang 
48642c57489SHong Zhang #undef __FUNCT__
48742c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
48842c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
48942c57489SHong Zhang {
49042c57489SHong Zhang   PetscErrorCode       ierr;
49142c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
492de0260b3SHong Zhang   Mat_MatMatMultMPI    *ap;
493776b82aeSLisandro Dalcin   PetscContainer       cont_merge,cont_ptap;
494de0260b3SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
49542c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
496de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
49742c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
498bf3909cdSBarry Smith   PetscInt             *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp;
49982412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
50082412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
50182412ba7SHong Zhang   MatScalar            *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth;
502d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5037adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
50442c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
50582412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
50642c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
50750f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
50842c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
50942c57489SHong Zhang   MPI_Status           *status;
51082412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
51182412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
512d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
51342c57489SHong Zhang 
51442c57489SHong Zhang   PetscFunctionBegin;
51542c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
51642c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
51742c57489SHong Zhang 
51842c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr);
51942c57489SHong Zhang   if (cont_merge) {
520776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr);
52142c57489SHong Zhang   } else {
52242c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
52342c57489SHong Zhang   }
5246c6e00e0SHong Zhang 
525f33d1a9aSHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr);
52642c57489SHong Zhang   if (cont_ptap) {
527776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr);
5286c6e00e0SHong Zhang     if (ap->reuse == MAT_INITIAL_MATRIX){
5296c6e00e0SHong Zhang       ap->reuse = MAT_REUSE_MATRIX;
5306c6e00e0SHong Zhang     } else { /* update numerical values of P_oth and P_loc */
5311d79065fSBarry Smith       ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
5326c6e00e0SHong Zhang       ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr);
5336c6e00e0SHong Zhang     }
53442c57489SHong Zhang   } else {
53542c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
53642c57489SHong Zhang   }
53742c57489SHong Zhang 
53842c57489SHong Zhang   /* get data from symbolic products */
539de0260b3SHong Zhang   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
540de0260b3SHong Zhang   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
54182412ba7SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc;
54242c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
54342c57489SHong Zhang 
544de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
545de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
546de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
547de0260b3SHong Zhang 
548de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5497a2fc3feSBarry Smith   owners = merge->rowmap->range;
550de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
551de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
55242c57489SHong Zhang 
553f4a743e1SHong Zhang   /* get data from symbolic A*P */
554de0260b3SHong Zhang   ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
555de0260b3SHong Zhang   ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
55642c57489SHong Zhang 
557f4a743e1SHong Zhang   /* compute numeric C_seq=P_loc^T*A_loc*P */
55882412ba7SHong Zhang   api = ap->abi; apj = ap->abj;
55942c57489SHong Zhang   for (i=0;i<am;i++) {
560f4a743e1SHong Zhang     /* form i-th sparse row of A*P */
56182412ba7SHong Zhang     apnz = api[i+1] - api[i];
56282412ba7SHong Zhang     apJ  = apj + api[i];
56342c57489SHong Zhang     /* diagonal portion of A */
56482412ba7SHong Zhang     anz  = adi[i+1] - adi[i];
56582412ba7SHong Zhang     for (j=0;j<anz;j++) {
56682412ba7SHong Zhang       row = *adj++;
56782412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
56882412ba7SHong Zhang       pj  = pj_loc + pi_loc[row];
56982412ba7SHong Zhang       pa  = pa_loc + pi_loc[row];
570f4a743e1SHong Zhang       nextp = 0;
57182412ba7SHong Zhang       for (k=0; nextp<pnz; k++) {
57282412ba7SHong Zhang         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
57382412ba7SHong Zhang           apa[k] += (*ada)*pa[nextp++];
57442c57489SHong Zhang         }
57542c57489SHong Zhang       }
576dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
57742c57489SHong Zhang       ada++;
57842c57489SHong Zhang     }
57942c57489SHong Zhang     /* off-diagonal portion of A */
58082412ba7SHong Zhang     anz  = aoi[i+1] - aoi[i];
58182412ba7SHong Zhang     for (j=0; j<anz; j++) {
58282412ba7SHong Zhang       row = *aoj++;
58382412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
58482412ba7SHong Zhang       pj  = pj_oth + pi_oth[row];
58582412ba7SHong Zhang       pa  = pa_oth + pi_oth[row];
586f4a743e1SHong Zhang       nextp = 0;
58782412ba7SHong Zhang       for (k=0; nextp<pnz; k++) {
58882412ba7SHong Zhang         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
58982412ba7SHong Zhang           apa[k] += (*aoa)*pa[nextp++];
59042c57489SHong Zhang         }
59142c57489SHong Zhang       }
592dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
59342c57489SHong Zhang       aoa++;
59442c57489SHong Zhang     }
59542c57489SHong Zhang 
59682412ba7SHong Zhang     /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
59782412ba7SHong Zhang     pnz = pi_loc[i+1] - pi_loc[i];
59882412ba7SHong Zhang     for (j=0; j<pnz; j++) {
59942c57489SHong Zhang       nextap = 0;
60082412ba7SHong Zhang       row    = *pJ++; /* global index */
60182412ba7SHong Zhang       if (row < pcstart || row >=pcend) { /* put the value into Co */
60282412ba7SHong Zhang         cj  = coj + coi[*poJ];
60382412ba7SHong Zhang         ca  = coa + coi[*poJ++];
60482412ba7SHong Zhang       } else {                            /* put the value into Cd */
60582412ba7SHong Zhang         cj   = bj + bi[*pdJ];
60682412ba7SHong Zhang         ca   = ba + bi[*pdJ++];
60742c57489SHong Zhang       }
60882412ba7SHong Zhang       for (k=0; nextap<apnz; k++) {
60982412ba7SHong Zhang         if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++];
61042c57489SHong Zhang       }
611dc0b31edSSatish Balay       ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
61282412ba7SHong Zhang       pA++;
613de0260b3SHong Zhang     }
614de0260b3SHong Zhang 
61542c57489SHong Zhang     /* zero the current row info for A*P */
61682412ba7SHong Zhang     ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
61742c57489SHong Zhang   }
61842c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
61942c57489SHong Zhang 
620de0260b3SHong Zhang   /* send and recv matrix values */
621de0260b3SHong Zhang   /*-----------------------------*/
62242c57489SHong Zhang   buf_ri = merge->buf_ri;
62342c57489SHong Zhang   buf_rj = merge->buf_rj;
62442c57489SHong Zhang   len_s  = merge->len_s;
625fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
62642c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
62742c57489SHong Zhang 
62842c57489SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr);
62942c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
63042c57489SHong Zhang     if (!len_s[proc]) continue;
631de0260b3SHong Zhang     i = merge->owners_co[proc];
632de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
63342c57489SHong Zhang     k++;
63442c57489SHong Zhang   }
63582412ba7SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
6360c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
6370c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
63842c57489SHong Zhang   ierr = PetscFree(status);CHKERRQ(ierr);
63942c57489SHong Zhang 
64042c57489SHong Zhang   ierr = PetscFree(s_waits);CHKERRQ(ierr);
64142c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
64282412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
64342c57489SHong Zhang 
64482412ba7SHong Zhang   /* insert local and received values into C */
64582412ba7SHong Zhang   /*-----------------------------------------*/
646687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
64742c57489SHong Zhang 
64842c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
64942c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
65042c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
65142c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
65282412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
65342c57489SHong Zhang   }
65442c57489SHong Zhang 
655de0260b3SHong Zhang   for (i=0; i<cm; i++) {
65682412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
65742c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
658de0260b3SHong Zhang     ba_i = ba + bi[i];
65982412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
66042c57489SHong Zhang     /* add received vals into ba */
66142c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
66242c57489SHong Zhang       /* i-th row */
66342c57489SHong Zhang       if (i == *nextrow[k]) {
66482412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
66582412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
66682412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
66782412ba7SHong Zhang         nextcj = 0;
66882412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
66982412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
67082412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
67142c57489SHong Zhang           }
67242c57489SHong Zhang         }
67382412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
67442c57489SHong Zhang       }
67542c57489SHong Zhang     }
67682412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
677dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
67842c57489SHong Zhang   }
679fd0ff01cSHong Zhang   ierr = MatSetBlockSize(C,1);CHKERRQ(ierr);
68042c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68142c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68242c57489SHong Zhang 
683de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
684*c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
68542c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
6860572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
68742c57489SHong Zhang   PetscFunctionReturn(0);
68842c57489SHong Zhang }
689