xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 38571add2b95ecffbc5f70117e513581cc74908d)
182412ba7SHong Zhang 
242c57489SHong Zhang /*
342c57489SHong Zhang   Defines projective product routines where A is a MPIAIJ matrix
442c57489SHong Zhang           C = P^T * A * P
542c57489SHong Zhang */
642c57489SHong Zhang 
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 <../src/mat/impls/aij/mpi/mpiaij.h>
10c6db04a5SJed Brown #include <petscbt.h>
1142c57489SHong Zhang 
1224ecddacSHong Zhang /* #define PTAP_PROFILE */
1324ecddacSHong Zhang 
1409573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
15cc6cb787SHong Zhang #undef __FUNCT__
16f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
17f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
18cc6cb787SHong Zhang {
19cc6cb787SHong Zhang   PetscErrorCode       ierr;
20f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
21f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
22cc6cb787SHong Zhang 
23cc6cb787SHong Zhang   PetscFunctionBegin;
24f8487c73SHong Zhang   if (ptap){
25c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
26b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
27f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
28a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
29a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
30c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
31c5af039cSHong Zhang     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
32c5af039cSHong Zhang     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
33414894bdSHong Zhang     if (ptap->apa){ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
34c8b0795cSMark F. Adams     if (merge) {
35cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
36cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
37cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
38cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
39cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
40c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
41cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
42c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
43cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4405b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4505b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4605b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
476bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
48dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
49f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
50bf0cc555SLisandro Dalcin     }
51f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
52c8b0795cSMark F. Adams   }
53cc6cb787SHong Zhang   PetscFunctionReturn(0);
54cc6cb787SHong Zhang }
55cc6cb787SHong Zhang 
56cc6cb787SHong Zhang #undef __FUNCT__
57cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
58f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
59f0c0a3a6SBarry Smith {
60cc6cb787SHong Zhang   PetscErrorCode       ierr;
6111a6856fSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
6211a6856fSHong Zhang   Mat_PtAPMPI          *ptap = a->ptap;
6311a6856fSHong Zhang   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
64cc6cb787SHong Zhang 
65cc6cb787SHong Zhang   PetscFunctionBegin;
66dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
6711a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
6811a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
69cc6cb787SHong Zhang   PetscFunctionReturn(0);
70cc6cb787SHong Zhang }
71cc6cb787SHong Zhang 
7242c57489SHong Zhang #undef __FUNCT__
73cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
74cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
7542c57489SHong Zhang {
7642c57489SHong Zhang   PetscErrorCode ierr;
7742c57489SHong Zhang 
7842c57489SHong Zhang   PetscFunctionBegin;
79cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
80cf3ca8ceSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
81cf3ca8ceSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
82cf3ca8ceSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
837a7894deSKris Buschelman   }
84cf3ca8ceSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
85cf3ca8ceSHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
86cf3ca8ceSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
8742c57489SHong Zhang   PetscFunctionReturn(0);
8842c57489SHong Zhang }
8942c57489SHong Zhang 
9042c57489SHong Zhang #undef __FUNCT__
9142c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9242c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9342c57489SHong Zhang {
9442c57489SHong Zhang   PetscErrorCode       ierr;
9577584fe6SHong Zhang   Mat                  Cmpi;
96f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
97a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
98f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
9942c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10042c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
101ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10277584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
103a914927fSHong Zhang   PetscInt             *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
104d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
10542c57489SHong Zhang   PetscBT              lnkbt;
1067adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
107a914927fSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
10842c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
10942c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
11024ecddacSHong Zhang   PetscInt             nzi,*pti,*ptj;
11142c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
112ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
113ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11442c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
11577584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
116d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
117a914927fSHong Zhang   PetscInt             rmax;
11824ecddacSHong Zhang #if defined(PTAP_PROFILE)
119*38571addSHong Zhang   PetscLogDouble       t0,t1,t2,t3,t4;
12024ecddacSHong Zhang #endif
12142c57489SHong Zhang 
12242c57489SHong Zhang   PetscFunctionBegin;
12324ecddacSHong Zhang #if defined(PTAP_PROFILE)
12424ecddacSHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
12524ecddacSHong Zhang #endif
12624ecddacSHong Zhang 
127c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
128c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
129c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
130c5af039cSHong Zhang   }
131c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
132c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
133c5af039cSHong Zhang   }
134c5af039cSHong Zhang 
13583445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13683445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13783445d95SHong Zhang 
138f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
139f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
140f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
1416c6e00e0SHong Zhang 
1426c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
143b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
144fe615159SHong Zhang 
1456c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
146a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1476c6e00e0SHong Zhang 
148a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
149a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1506c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1516c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
15224ecddacSHong Zhang #if defined(PTAP_PROFILE)
15324ecddacSHong Zhang   ierr = PetscGetTime(&t1);CHKERRQ(ierr);
15424ecddacSHong Zhang #endif
15542c57489SHong Zhang 
156fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
157fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
15877584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
15982412ba7SHong Zhang   api[0] = 0;
16083445d95SHong Zhang 
16183445d95SHong Zhang   /* create and initialize a linked list */
162a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
16383445d95SHong Zhang 
164a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
165d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
166f4a743e1SHong Zhang   current_space = free_space;
167f4a743e1SHong Zhang 
168f4a743e1SHong Zhang   for (i=0; i<am; i++) {
169f4a743e1SHong Zhang     /* diagonal portion of A */
170ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
17177584fe6SHong Zhang     aj  = ad->j + adi[i];
172ded4f561SHong Zhang     for (j=0; j<nzi; j++){
17377584fe6SHong Zhang       row  = aj[j];
17482412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
175ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
17683445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
177a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
178f4a743e1SHong Zhang     }
179f4a743e1SHong Zhang     /* off-diagonal portion of A */
180ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
18177584fe6SHong Zhang     aj  = ao->j + aoi[i];
182ded4f561SHong Zhang     for (j=0; j<nzi; j++){
18377584fe6SHong Zhang       row = aj[j];
18482412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
185ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
186a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
187f4a743e1SHong Zhang     }
188a914927fSHong Zhang     apnz = lnk[0];
18982412ba7SHong Zhang     api[i+1] = api[i] + apnz;
19077584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
191f4a743e1SHong Zhang 
19283445d95SHong Zhang     /* if free space is not available, double the total space in the list */
19382412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1942ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
195f2b054eeSHong Zhang       nspacedouble++;
196f4a743e1SHong Zhang     }
197f4a743e1SHong Zhang 
198f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
199a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
20082412ba7SHong Zhang     current_space->array           += apnz;
20182412ba7SHong Zhang     current_space->local_used      += apnz;
20282412ba7SHong Zhang     current_space->local_remaining -= apnz;
203f4a743e1SHong Zhang   }
204a914927fSHong Zhang 
20582412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
206f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
20777584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
20877584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
209118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
210d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
211f4a743e1SHong Zhang 
21224ecddacSHong Zhang #if defined(PTAP_PROFILE)
21324ecddacSHong Zhang   ierr = PetscGetTime(&t2);CHKERRQ(ierr);
21424ecddacSHong Zhang #endif
21524ecddacSHong 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 
226d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
227d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
228a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
22942c57489SHong Zhang   current_space = free_space;
23042c57489SHong Zhang 
231de0260b3SHong Zhang   for (i=0; i<pon; i++) {
23282412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
23377584fe6SHong Zhang     ptJ = potj + poti[i];
23477584fe6SHong Zhang     for (j=0; j<pnz; j++){
23577584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
23682412ba7SHong Zhang       apnz = api[row+1] - api[row];
237ded4f561SHong Zhang       Jptr = apj + api[row];
23883445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
239a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
24042c57489SHong Zhang     }
241a914927fSHong Zhang     nnz = lnk[0];
24242c57489SHong Zhang 
24383445d95SHong Zhang     /* If free space is not available, double the total space in the list */
244ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2454238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
246d16ca5f0SHong Zhang       nspacedouble++;
24742c57489SHong Zhang     }
24842c57489SHong Zhang 
24942c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
250a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
251ded4f561SHong Zhang     current_space->array           += nnz;
252ded4f561SHong Zhang     current_space->local_used      += nnz;
253ded4f561SHong Zhang     current_space->local_remaining -= nnz;
254ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
25542c57489SHong Zhang   }
256de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
257a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
258118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
259d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
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 */
27324ecddacSHong Zhang   ierr = PetscMalloc2(size,PetscMPIInt,&len_si,size,MPI_Status,&sstatus);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 */
305529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
306529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
307ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
30842c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
30942c57489SHong Zhang     if (!len_s[proc]) continue;
310de0260b3SHong Zhang     i = owners_co[proc];
311529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
31242c57489SHong Zhang     k++;
31342c57489SHong Zhang   }
31442c57489SHong Zhang 
315ded4f561SHong Zhang   /* receives and sends of coj are complete */
31677584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
317c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
318ded4f561SHong Zhang   }
319ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3200c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
32182412ba7SHong Zhang 
322ded4f561SHong Zhang   /* send and recv coi */
323ded4f561SHong Zhang   /*-------------------*/
324529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
325529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
32642c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
32742c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
32842c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
32942c57489SHong Zhang     if (!len_s[proc]) continue;
33042c57489SHong Zhang     /* form outgoing message for i-structure:
33142c57489SHong Zhang          buf_si[0]:                 nrows to be sent
33242c57489SHong Zhang                [1:nrows]:           row index (global)
33342c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
33442c57489SHong Zhang     */
33542c57489SHong Zhang     /*-------------------------------------------*/
33642c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
33742c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
33842c57489SHong Zhang     buf_si[0]   = nrows;
33942c57489SHong Zhang     buf_si_i[0] = 0;
34042c57489SHong Zhang     nrows = 0;
341de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
342ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
343ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
344de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
34542c57489SHong Zhang       nrows++;
34642c57489SHong Zhang     }
347529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
34842c57489SHong Zhang     k++;
34942c57489SHong Zhang     buf_si += len_si[proc];
35042c57489SHong Zhang   }
351ded4f561SHong Zhang   i = merge->nrecv;
352ded4f561SHong Zhang   while (i--) {
353c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
354ded4f561SHong Zhang   }
355ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3560c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
357a914927fSHong Zhang 
35824ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
35942c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
360ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
36142c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
36242c57489SHong Zhang 
363*38571addSHong Zhang #if defined(PTAP_PROFILE)
364*38571addSHong Zhang   ierr = PetscGetTime(&t3);CHKERRQ(ierr);
365*38571addSHong Zhang #endif
366*38571addSHong Zhang 
367de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
368de0260b3SHong Zhang   /*------------------------------------------*/
369ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
370ded4f561SHong Zhang 
37124ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
37224ecddacSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&pti);CHKERRQ(ierr);
37324ecddacSHong Zhang   pti[0] = 0;
37442c57489SHong Zhang 
375d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
376d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
377a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
37842c57489SHong Zhang   current_space = free_space;
37942c57489SHong Zhang 
380687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
38142c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
38242c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
38342c57489SHong Zhang     nrows       = *buf_ri_k[k];
38442c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
38542c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
38642c57489SHong Zhang   }
38742c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
38808cb4532SHong Zhang   rmax = 0;
38942c57489SHong Zhang   for (i=0; i<pn; i++) {
390ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
391ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
39277584fe6SHong Zhang     ptJ = pdtj + pdti[i];
39377584fe6SHong Zhang     for (j=0; j<pnz; j++){
39477584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
395ded4f561SHong Zhang       apnz = api[row+1] - api[row];
396ded4f561SHong Zhang       Jptr = apj + api[row];
397ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
398a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
399ded4f561SHong Zhang     }
400a914927fSHong Zhang 
40142c57489SHong Zhang     /* add received col data into lnk */
40242c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
40342c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
404ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
405ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
406a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
40742c57489SHong Zhang         nextrow[k]++; nextci[k]++;
40842c57489SHong Zhang       }
40942c57489SHong Zhang     }
410a914927fSHong Zhang     nnz = lnk[0];
41142c57489SHong Zhang 
41242c57489SHong Zhang     /* if free space is not available, make more free space */
413ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4144238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
415d16ca5f0SHong Zhang       nspacedouble++;
41642c57489SHong Zhang     }
41742c57489SHong Zhang     /* copy data into free space, then initialize lnk */
418a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
419ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
420ded4f561SHong Zhang     current_space->array           += nnz;
421ded4f561SHong Zhang     current_space->local_used      += nnz;
422ded4f561SHong Zhang     current_space->local_remaining -= nnz;
42324ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
42408cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
42542c57489SHong Zhang   }
426ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4270572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
42842c57489SHong Zhang 
42924ecddacSHong Zhang   ierr = PetscMalloc((pti[pn]+1)*sizeof(PetscInt),&ptj);CHKERRQ(ierr);
43024ecddacSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
43124ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
432d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
43342c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
43442c57489SHong Zhang 
435a914927fSHong Zhang   /* create symbolic parallel matrix Cmpi */
436a914927fSHong Zhang   /*--------------------------------------*/
43777584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
43877584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
439a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
44077584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
44177584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
44242c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
443a2f3521dSMark F. Adams 
44424ecddacSHong Zhang   merge->bi            = pti;
44524ecddacSHong Zhang   merge->bj            = ptj;
446de0260b3SHong Zhang   merge->coi           = coi;
447de0260b3SHong Zhang   merge->coj           = coj;
44842c57489SHong Zhang   merge->buf_ri        = buf_ri;
44942c57489SHong Zhang   merge->buf_rj        = buf_rj;
450de0260b3SHong Zhang   merge->owners_co     = owners_co;
45177584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
45277584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
453cc6cb787SHong Zhang 
45477584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
455a914927fSHong Zhang   Cmpi->assembled      = PETSC_FALSE;
45677584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
45777584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
45842c57489SHong Zhang 
45977584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
46077584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
461f8487c73SHong Zhang   c->ptap        = ptap;
46277584fe6SHong Zhang   ptap->api      = api;
46377584fe6SHong Zhang   ptap->apj      = apj;
464f8487c73SHong Zhang   ptap->merge    = merge;
465d6ab1dc0SHong Zhang   ptap->rmax     = ap_rmax;
46677584fe6SHong Zhang   *C             = Cmpi;
467414894bdSHong Zhang 
468414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
469414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
470414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
471414894bdSHong Zhang   /* set default scalable */
472414894bdSHong Zhang   ptap->scalable = PETSC_TRUE;
473414894bdSHong Zhang   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,PETSC_NULL);CHKERRQ(ierr);
474414894bdSHong Zhang   if (!ptap->scalable){  /* Do dense axpy */
475414894bdSHong Zhang     ierr = PetscMalloc(pN*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
476*38571addSHong Zhang     ierr = PetscMemzero(ptap->apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
477414894bdSHong Zhang   } else {
478414894bdSHong Zhang     ierr = PetscMalloc((ap_rmax+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
479*38571addSHong Zhang     ierr = PetscMemzero(ptap->apa,ap_rmax*sizeof(PetscScalar));CHKERRQ(ierr);
480414894bdSHong Zhang   }
481414894bdSHong Zhang 
48224ecddacSHong Zhang #if defined(PTAP_PROFILE)
483*38571addSHong Zhang   ierr = PetscGetTime(&t4);CHKERRQ(ierr);
484*38571addSHong Zhang   if (rank==1) PetscPrintf(MPI_COMM_SELF,"  [%d] PtAPSymbolic %g/P + %g/AP + %g/comm + %g/PtAP = %g\n",rank,t1-t0,t2-t1,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr);
48524ecddacSHong Zhang #endif
48624ecddacSHong Zhang 
487f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
48824ecddacSHong Zhang   if (pti[pn] != 0) {
48977584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
49077584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
491f2b054eeSHong Zhang   } else {
49277584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
493f2b054eeSHong Zhang   }
494f2b054eeSHong Zhang #endif
49542c57489SHong Zhang   PetscFunctionReturn(0);
49642c57489SHong Zhang }
49742c57489SHong Zhang 
49842c57489SHong Zhang #undef __FUNCT__
49942c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
50042c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
50142c57489SHong Zhang {
50242c57489SHong Zhang   PetscErrorCode       ierr;
50342c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
504f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
50542c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
506de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
50742c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
508f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
5091d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
51082412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
51182412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
512e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
513d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5147adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
51542c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
51682412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
51742c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
51850f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
51942c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
52042c57489SHong Zhang   MPI_Status           *status;
52182412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
52282412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
523d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
5246ce94e8fSHong Zhang   PetscBool            scalable;
525*38571addSHong Zhang #if defined(PTAP_PROFILE)
526*38571addSHong Zhang   PetscLogDouble       t0,t1,t2,t3,t4;
527*38571addSHong Zhang #endif
52842c57489SHong Zhang 
52942c57489SHong Zhang   PetscFunctionBegin;
530*38571addSHong Zhang #if defined(PTAP_PROFILE)
531*38571addSHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
532*38571addSHong Zhang #endif
53342c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
53442c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
53542c57489SHong Zhang 
536f8487c73SHong Zhang   ptap = c->ptap;
5371c4805f8SJed Brown   if (!ptap) SETERRQ(((PetscObject)C)->comm,PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
538f8487c73SHong Zhang   merge    = ptap->merge;
539414894bdSHong Zhang   apa      = ptap->apa;
5406ce94e8fSHong Zhang   scalable = ptap->scalable;
5416c6e00e0SHong Zhang 
542a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
543f9473cf7SHong Zhang   /*--------------------------------------------------*/
544f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
545f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5466c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
547b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
548a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5496c6e00e0SHong Zhang   }
550*38571addSHong Zhang #if defined(PTAP_PROFILE)
551*38571addSHong Zhang   ierr = PetscGetTime(&t1);CHKERRQ(ierr);
552*38571addSHong Zhang #endif
553f8487c73SHong Zhang 
554f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
555f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
55642c57489SHong Zhang   /* get data from symbolic products */
557a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
558a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
559a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
56042c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
56142c57489SHong Zhang 
562de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
563de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
564de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
565de0260b3SHong Zhang 
566de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5677a2fc3feSBarry Smith   owners = merge->rowmap->range;
568de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
569de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
57042c57489SHong Zhang 
571a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
572d9824c15SHong Zhang 
573*38571addSHong Zhang   if (!scalable){ /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but fast */
574*38571addSHong Zhang     /*-----------------------------------------------------------------------------------------------------*/
575a9d06231SHong Zhang     for (i=0; i<am; i++) {
576a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
577a9d06231SHong Zhang       /*------------------------------------------------------------*/
578a9d06231SHong Zhang       apJ = apj + api[i];
579a9d06231SHong Zhang 
580a9d06231SHong Zhang       /* diagonal portion of A */
581a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
582a9d06231SHong Zhang       adj = ad->j + adi[i];
583a9d06231SHong Zhang       ada = ad->a + adi[i];
584a9d06231SHong Zhang       for (j=0; j<anz; j++) {
585a9d06231SHong Zhang         row = adj[j];
586a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
587a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
588a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
589a9d06231SHong Zhang 
590a9d06231SHong Zhang         /* perform dense axpy */
591e900a757SHong Zhang         valtmp = ada[j];
592a9d06231SHong Zhang         for (k=0; k<pnz; k++){
593e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
594a9d06231SHong Zhang         }
595a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
596a9d06231SHong Zhang       }
597a9d06231SHong Zhang 
598a9d06231SHong Zhang       /* off-diagonal portion of A */
599a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
600a9d06231SHong Zhang       aoj = ao->j + aoi[i];
601a9d06231SHong Zhang       aoa = ao->a + aoi[i];
602a9d06231SHong Zhang       for (j=0; j<anz; j++) {
603a9d06231SHong Zhang         row = aoj[j];
604a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
605a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
606a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
607a9d06231SHong Zhang 
608a9d06231SHong Zhang         /* perform dense axpy */
609e900a757SHong Zhang         valtmp = aoa[j];
610a9d06231SHong Zhang         for (k=0; k<pnz; k++){
611e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
612a9d06231SHong Zhang         }
613a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
614a9d06231SHong Zhang       }
615a9d06231SHong Zhang 
616a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
617a9d06231SHong Zhang       /*--------------------------------------------------------------*/
618a9d06231SHong Zhang       apnz = api[i+1] - api[i];
619a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
620a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
621e900a757SHong Zhang       poJ = po->j + po->i[i];
622a9d06231SHong Zhang       pA  = po->a + po->i[i];
623a9d06231SHong Zhang       for (j=0; j<pnz; j++){
624e900a757SHong Zhang         row = poJ[j];
625e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
626e900a757SHong Zhang         cj  = coj + coi[row];
627e900a757SHong Zhang         ca  = coa + coi[row];
628a9d06231SHong Zhang         /* perform dense axpy */
629e900a757SHong Zhang         valtmp = pA[j];
630a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
631e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
632a9d06231SHong Zhang         }
633a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
634a9d06231SHong Zhang       }
635a9d06231SHong Zhang 
636a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
637a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
638e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
639a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
640a9d06231SHong Zhang       for (j=0; j<pnz; j++){
641e900a757SHong Zhang         row = pdJ[j];
642e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
643e900a757SHong Zhang         cj  = bj + bi[row];
644e900a757SHong Zhang         ca  = ba + bi[row];
645a9d06231SHong Zhang         /* perform dense axpy */
646e900a757SHong Zhang         valtmp = pA[j];
647a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
648e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
649a9d06231SHong Zhang         }
650a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
651a9d06231SHong Zhang       }
652a9d06231SHong Zhang 
653a9d06231SHong Zhang       /* zero the current row of A*P */
654a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
655a9d06231SHong Zhang     }
656*38571addSHong Zhang   } else {/* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
657*38571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
658a9d06231SHong Zhang     pA=pa_loc;
65942c57489SHong Zhang     for (i=0; i<am; i++) {
660f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
66182412ba7SHong Zhang       apnz = api[i+1] - api[i];
66282412ba7SHong Zhang       apJ  = apj + api[i];
66342c57489SHong Zhang       /* diagonal portion of A */
66482412ba7SHong Zhang       anz = adi[i+1] - adi[i];
6651d633602SHong Zhang       adj = ad->j + adi[i];
6661d633602SHong Zhang       ada = ad->a + adi[i];
66782412ba7SHong Zhang       for (j=0; j<anz; j++) {
6681d633602SHong Zhang         row = adj[j];
66982412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
67082412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
67182412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
672e900a757SHong Zhang         valtmp = ada[j];
673f4a743e1SHong Zhang         nextp  = 0;
67482412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
67582412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
676e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
67742c57489SHong Zhang           }
67842c57489SHong Zhang         }
679dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
68042c57489SHong Zhang       }
68142c57489SHong Zhang       /* off-diagonal portion of A */
68282412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
6831d633602SHong Zhang       aoj = ao->j + aoi[i];
6841d633602SHong Zhang       aoa = ao->a + aoi[i];
68582412ba7SHong Zhang       for (j=0; j<anz; j++) {
6861d633602SHong Zhang         row = aoj[j];
68782412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
68882412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
68982412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
690e900a757SHong Zhang         valtmp = aoa[j];
691f4a743e1SHong Zhang         nextp  = 0;
69282412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
69382412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
694e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
69542c57489SHong Zhang           }
69642c57489SHong Zhang         }
697dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
69842c57489SHong Zhang       }
69942c57489SHong Zhang 
700a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
701a9d06231SHong Zhang       /*--------------------------------------------------------------*/
70282412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
703f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
70482412ba7SHong Zhang       for (j=0; j<pnz; j++) {
70542c57489SHong Zhang         nextap = 0;
706f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
70782412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
708e900a757SHong Zhang           row = *poJ;
709e900a757SHong Zhang           cj  = coj + coi[row];
710e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
71182412ba7SHong Zhang         } else {                            /* put the value into Cd */
712e900a757SHong Zhang           row  = *pdJ;
713e900a757SHong Zhang           cj   = bj + bi[row];
714e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
71542c57489SHong Zhang         }
716e900a757SHong Zhang         valtmp = pA[j];
71782412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
718e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
71942c57489SHong Zhang         }
720dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
721de0260b3SHong Zhang       }
7220496f32cSHong Zhang       pA += pnz;
72342c57489SHong Zhang       /* zero the current row info for A*P */
72482412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
72542c57489SHong Zhang     }
726d9824c15SHong Zhang   }
727*38571addSHong Zhang #if defined(PTAP_PROFILE)
728*38571addSHong Zhang   ierr = PetscGetTime(&t2);CHKERRQ(ierr);
729*38571addSHong Zhang #endif
73042c57489SHong Zhang 
731a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
732a9d06231SHong Zhang   /*------------------------------------*/
73342c57489SHong Zhang   buf_ri = merge->buf_ri;
73442c57489SHong Zhang   buf_rj = merge->buf_rj;
73542c57489SHong Zhang   len_s  = merge->len_s;
736fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
73742c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
73842c57489SHong Zhang 
739a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
74042c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
74142c57489SHong Zhang     if (!len_s[proc]) continue;
742de0260b3SHong Zhang     i = merge->owners_co[proc];
743de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
74442c57489SHong Zhang     k++;
74542c57489SHong Zhang   }
7460c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
7470c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
74842c57489SHong Zhang 
749a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
75042c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
75182412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
752*38571addSHong Zhang #if defined(PTAP_PROFILE)
753*38571addSHong Zhang   ierr = PetscGetTime(&t3);CHKERRQ(ierr);
754*38571addSHong Zhang #endif
75542c57489SHong Zhang 
756a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
757a9d06231SHong Zhang   /*------------------------------------------------------*/
758687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
75942c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
76042c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
76142c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
76242c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
76382412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
76442c57489SHong Zhang   }
76542c57489SHong Zhang 
766de0260b3SHong Zhang   for (i=0; i<cm; i++) {
76782412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
76842c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
769de0260b3SHong Zhang     ba_i = ba + bi[i];
77082412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
77142c57489SHong Zhang     /* add received vals into ba */
77242c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
77342c57489SHong Zhang       /* i-th row */
77442c57489SHong Zhang       if (i == *nextrow[k]) {
77582412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
77682412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
77782412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
77882412ba7SHong Zhang         nextcj = 0;
77982412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
78082412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
78182412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
78242c57489SHong Zhang           }
78342c57489SHong Zhang         }
78482412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
785c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
78642c57489SHong Zhang       }
78742c57489SHong Zhang     }
78882412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
78942c57489SHong Zhang   }
79042c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79142c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79242c57489SHong Zhang 
793de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
794c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
79542c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
7960572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
797*38571addSHong Zhang #if defined(PTAP_PROFILE)
798*38571addSHong Zhang   ierr = PetscGetTime(&t4);CHKERRQ(ierr);
799*38571addSHong Zhang   if (rank==1) PetscPrintf(MPI_COMM_SELF,"  [%d] PtAPNum %g/P + %g/PtAP + %g/comm + %g/Cloc = %g\n\n",rank,t1-t0,t2-t1,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr);
800*38571addSHong Zhang #endif
801*38571addSHong Zhang 
802*38571addSHong Zhang 
803d9824c15SHong Zhang #if defined(PETSC_USE_INFO)
804d9824c15SHong Zhang   if (scalable){
805d9824c15SHong Zhang     ierr = PetscInfo(C,"Use scalable sparse axpy\n");CHKERRQ(ierr);
806d9824c15SHong Zhang   } else {
807d9824c15SHong Zhang     ierr = PetscInfo(C,"Use non-scalable dense axpy\n");CHKERRQ(ierr);
808d9824c15SHong Zhang   }
809d9824c15SHong Zhang #endif
81042c57489SHong Zhang   PetscFunctionReturn(0);
81142c57489SHong Zhang }
812