xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 6b911d1655e60e81bc8952dcf2b0a2eb7516084c)
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>
118563dfccSBarry Smith #include <petsctime.h>
1242c57489SHong Zhang 
1324ecddacSHong Zhang /* #define PTAP_PROFILE */
1424ecddacSHong Zhang 
1509573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
16cc6cb787SHong Zhang #undef __FUNCT__
17f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
18f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
19cc6cb787SHong Zhang {
20cc6cb787SHong Zhang   PetscErrorCode ierr;
21f8487c73SHong Zhang   Mat_MPIAIJ     *a   =(Mat_MPIAIJ*)A->data;
22f8487c73SHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
23cc6cb787SHong Zhang 
24cc6cb787SHong Zhang   PetscFunctionBegin;
25f8487c73SHong Zhang   if (ptap) {
26c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI *merge=ptap->merge;
27b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
28f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
29a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
30a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
31c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
32c5af039cSHong Zhang     if (ptap->api) {ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
33c5af039cSHong Zhang     if (ptap->apj) {ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
34414894bdSHong Zhang     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
35c8b0795cSMark F. Adams     if (merge) {
36cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
37cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
38cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
39cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
40cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
41c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
42cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
43c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
44cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4505b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4605b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4705b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
486bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
49dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
50f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
51bf0cc555SLisandro Dalcin     }
52f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
53c8b0795cSMark F. Adams   }
54cc6cb787SHong Zhang   PetscFunctionReturn(0);
55cc6cb787SHong Zhang }
56cc6cb787SHong Zhang 
57cc6cb787SHong Zhang #undef __FUNCT__
58cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
59f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
60f0c0a3a6SBarry Smith {
61cc6cb787SHong Zhang   PetscErrorCode      ierr;
6211a6856fSHong Zhang   Mat_MPIAIJ          *a     = (Mat_MPIAIJ*)A->data;
6311a6856fSHong Zhang   Mat_PtAPMPI         *ptap  = a->ptap;
6411a6856fSHong Zhang   Mat_Merge_SeqsToMPI *merge = ptap->merge;
65cc6cb787SHong Zhang 
66cc6cb787SHong Zhang   PetscFunctionBegin;
67dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
682205254eSKarl Rupp 
6911a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
7011a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
71cc6cb787SHong Zhang   PetscFunctionReturn(0);
72cc6cb787SHong Zhang }
73cc6cb787SHong Zhang 
7442c57489SHong Zhang #undef __FUNCT__
75cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
76cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
7742c57489SHong Zhang {
7842c57489SHong Zhang   PetscErrorCode ierr;
7942c57489SHong Zhang 
8042c57489SHong Zhang   PetscFunctionBegin;
81cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
823ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
83cf3ca8ceSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
843ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
857a7894deSKris Buschelman   }
863ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
87cf3ca8ceSHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
883ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
8942c57489SHong Zhang   PetscFunctionReturn(0);
9042c57489SHong Zhang }
9142c57489SHong Zhang 
9242c57489SHong Zhang #undef __FUNCT__
9342c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9442c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9542c57489SHong Zhang {
9642c57489SHong Zhang   PetscErrorCode      ierr;
9777584fe6SHong Zhang   Mat                 Cmpi;
98f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
990298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
100f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10142c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10242c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
103ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10477584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
105a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
106d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
10742c57489SHong Zhang   PetscBT             lnkbt;
108ce94432eSBarry Smith   MPI_Comm            comm;
109a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
11042c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
11142c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
11224ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
11342c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
114ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
115ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
11642c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
11777584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
118d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
119a914927fSHong Zhang   PetscInt            rmax;
12042c57489SHong Zhang 
12142c57489SHong Zhang   PetscFunctionBegin;
122ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
12324ecddacSHong Zhang 
124c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
125c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
126c5af039cSHong 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);
127c5af039cSHong Zhang   }
128c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
129c5af039cSHong 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);
130c5af039cSHong Zhang   }
131c5af039cSHong Zhang 
13283445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13383445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13483445d95SHong Zhang 
135f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
136b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
137b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
1389f4155fbSHong Zhang   ptap->merge = merge;
139f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
1406c6e00e0SHong Zhang 
1416c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
142b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
143fe615159SHong Zhang 
1446c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
145a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1466c6e00e0SHong Zhang 
147a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
148a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1496c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1506c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
15142c57489SHong Zhang 
152*6b911d16SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
153fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
154854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
15582412ba7SHong Zhang   api[0] = 0;
15683445d95SHong Zhang 
15783445d95SHong Zhang   /* create and initialize a linked list */
158a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
15983445d95SHong Zhang 
160a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
161d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
1622205254eSKarl Rupp 
163f4a743e1SHong Zhang   current_space = free_space;
164f4a743e1SHong Zhang 
165f4a743e1SHong Zhang   for (i=0; i<am; i++) {
166f4a743e1SHong Zhang     /* diagonal portion of A */
167ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
16877584fe6SHong Zhang     aj  = ad->j + adi[i];
169ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
17077584fe6SHong Zhang       row  = aj[j];
17182412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
172ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
17383445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
174a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
175f4a743e1SHong Zhang     }
176f4a743e1SHong Zhang     /* off-diagonal portion of A */
177ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
17877584fe6SHong Zhang     aj  = ao->j + aoi[i];
179ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
18077584fe6SHong Zhang       row  = aj[j];
18182412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
182ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
183a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
184f4a743e1SHong Zhang     }
185a914927fSHong Zhang     apnz     = lnk[0];
18682412ba7SHong Zhang     api[i+1] = api[i] + apnz;
18777584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
188f4a743e1SHong Zhang 
18983445d95SHong Zhang     /* if free space is not available, double the total space in the list */
19082412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1912ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
192f2b054eeSHong Zhang       nspacedouble++;
193f4a743e1SHong Zhang     }
194f4a743e1SHong Zhang 
195f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
196a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1972205254eSKarl Rupp 
19882412ba7SHong Zhang     current_space->array           += apnz;
19982412ba7SHong Zhang     current_space->local_used      += apnz;
20082412ba7SHong Zhang     current_space->local_remaining -= apnz;
201f4a743e1SHong Zhang   }
202a914927fSHong Zhang 
20382412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
204f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
205854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
20677584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
207118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
208d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
209f4a743e1SHong Zhang 
210*6b911d16SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others */
211*6b911d16SHong Zhang   /*--------------------------------------------------------*/
212de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
21342c57489SHong Zhang 
214ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
215d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
21683445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
217854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
218de0260b3SHong Zhang   coi[0] = 0;
21942c57489SHong Zhang 
220d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
221d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
22222d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
22342c57489SHong Zhang   current_space = free_space;
22442c57489SHong Zhang 
225de0260b3SHong Zhang   for (i=0; i<pon; i++) {
22682412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
22777584fe6SHong Zhang     ptJ = potj + poti[i];
22877584fe6SHong Zhang     for (j=0; j<pnz; j++) {
22977584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
23082412ba7SHong Zhang       apnz = api[row+1] - api[row];
231ded4f561SHong Zhang       Jptr = apj + api[row];
23283445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
233a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
23442c57489SHong Zhang     }
235a914927fSHong Zhang     nnz = lnk[0];
23642c57489SHong Zhang 
23783445d95SHong Zhang     /* If free space is not available, double the total space in the list */
238ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2394238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
240d16ca5f0SHong Zhang       nspacedouble++;
24142c57489SHong Zhang     }
24242c57489SHong Zhang 
24342c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
244a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
2452205254eSKarl Rupp 
246ded4f561SHong Zhang     current_space->array           += nnz;
247ded4f561SHong Zhang     current_space->local_used      += nnz;
248ded4f561SHong Zhang     current_space->local_remaining -= nnz;
2492205254eSKarl Rupp 
250ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
25142c57489SHong Zhang   }
252*6b911d16SHong Zhang   printf("[%d] coi[%d] = %d\n",rank,pon,coi[pon]); //coi[pon] can be 0!!!
253*6b911d16SHong Zhang 
254*6b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
255a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
256118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
257d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
258de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
25942c57489SHong Zhang 
260*6b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
261*6b911d16SHong Zhang   /*--------------------------------------------------*/
262*6b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
263*6b911d16SHong Zhang   len_s        = merge->len_s;
264*6b911d16SHong Zhang   merge->nsend = 0;
265*6b911d16SHong Zhang 
266*6b911d16SHong Zhang 
26742c57489SHong Zhang   /* determine row ownership */
26826283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2697a2fc3feSBarry Smith   merge->rowmap->n  = pn;
2707a2fc3feSBarry Smith   merge->rowmap->bs = 1;
2712205254eSKarl Rupp 
27226283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2737a2fc3feSBarry Smith   owners = merge->rowmap->range;
27442c57489SHong Zhang 
27542c57489SHong Zhang   /* determine the number of messages to send, their lengths */
276dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
27783445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
278854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
279de0260b3SHong Zhang 
28083445d95SHong Zhang   proc = 0;
281de0260b3SHong Zhang   for (i=0; i<pon; i++) {
282de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
283ee6e4b5bSHong Zhang     len_si[proc]++;  /* num of rows in Co(=Pt*AP) to be sent to [proc] -- could be empty row!!! */
284ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
28583445d95SHong Zhang   }
286de0260b3SHong Zhang 
287de0260b3SHong Zhang   len          = 0; /* max length of buf_si[] */
288de0260b3SHong Zhang   owners_co[0] = 0;
28942c57489SHong Zhang   for (proc=0; proc<size; proc++) {
290de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
291ee6e4b5bSHong Zhang     if (len_s[proc]) {
29242c57489SHong Zhang       merge->nsend++;
29383445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
29442c57489SHong Zhang       len         += len_si[proc];
29542c57489SHong Zhang     }
29642c57489SHong Zhang   }
29742c57489SHong Zhang 
298ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
2990298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
300*6b911d16SHong Zhang   printf("[%d] nsend %d, nrecv %d\n",rank,merge->nsend,merge->nrecv);
30142c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
30242c57489SHong Zhang 
303ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
304529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
305529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
306854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
30742c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
30842c57489SHong Zhang     if (!len_s[proc]) continue;
309de0260b3SHong Zhang     i    = owners_co[proc];
310529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
31142c57489SHong Zhang     k++;
31242c57489SHong Zhang   }
31342c57489SHong Zhang 
314ded4f561SHong Zhang   /* receives and sends of coj are complete */
31577584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
316c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
317ded4f561SHong Zhang   }
318ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3190c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
32082412ba7SHong Zhang 
321*6b911d16SHong Zhang   /* (4) send and recv coi */
322*6b911d16SHong Zhang   /*-----------------------*/
323529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
324529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
325854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
32642c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
32742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
32842c57489SHong Zhang     if (!len_s[proc]) continue;
32942c57489SHong Zhang     /* form outgoing message for i-structure:
33042c57489SHong Zhang          buf_si[0]:                 nrows to be sent
33142c57489SHong Zhang                [1:nrows]:           row index (global)
33242c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
33342c57489SHong Zhang     */
33442c57489SHong Zhang     /*-------------------------------------------*/
33542c57489SHong Zhang     nrows       = len_si[proc]/2 - 1;
33642c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
33742c57489SHong Zhang     buf_si[0]   = nrows;
33842c57489SHong Zhang     buf_si_i[0] = 0;
33942c57489SHong Zhang     nrows       = 0;
340de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
341ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
3422205254eSKarl Rupp 
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*6b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
364*6b911d16SHong Zhang   /*----------------------------------------------*/
365ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
366ded4f561SHong Zhang 
36724ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
368854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
36924ecddacSHong Zhang   pti[0] = 0;
37042c57489SHong Zhang 
371d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
372d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
37322d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
37442c57489SHong Zhang   current_space = free_space;
37542c57489SHong Zhang 
376dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
37742c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
37842c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
37942c57489SHong Zhang     nrows       = *buf_ri_k[k];
38042c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
38142c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
38242c57489SHong Zhang   }
38342c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
38408cb4532SHong Zhang   rmax = 0;
38542c57489SHong Zhang   for (i=0; i<pn; i++) {
386ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
387ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
38877584fe6SHong Zhang     ptJ = pdtj + pdti[i];
38977584fe6SHong Zhang     for (j=0; j<pnz; j++) {
39077584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
391ded4f561SHong Zhang       apnz = api[row+1] - api[row];
392ded4f561SHong Zhang       Jptr = apj + api[row];
393ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
394a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
395ded4f561SHong Zhang     }
396a914927fSHong Zhang 
39742c57489SHong Zhang     /* add received col data into lnk */
39842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
39942c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
400ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
401ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
402a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
40342c57489SHong Zhang         nextrow[k]++; nextci[k]++;
40442c57489SHong Zhang       }
40542c57489SHong Zhang     }
406a914927fSHong Zhang     nnz = lnk[0];
40742c57489SHong Zhang 
40842c57489SHong Zhang     /* if free space is not available, make more free space */
409ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4104238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
411d16ca5f0SHong Zhang       nspacedouble++;
41242c57489SHong Zhang     }
41342c57489SHong Zhang     /* copy data into free space, then initialize lnk */
414a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
415ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
4162205254eSKarl Rupp 
417ded4f561SHong Zhang     current_space->array           += nnz;
418ded4f561SHong Zhang     current_space->local_used      += nnz;
419ded4f561SHong Zhang     current_space->local_remaining -= nnz;
4202205254eSKarl Rupp 
42124ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
42208cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
42342c57489SHong Zhang   }
424ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4250572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
42642c57489SHong Zhang 
427854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
42824ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
42924ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
430d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
43142c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
43242c57489SHong Zhang 
433*6b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
434*6b911d16SHong Zhang   /*------------------------------------------*/
43577584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
43677584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
43733d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
43877584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
43977584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
44042c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
441a2f3521dSMark F. Adams 
442b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
443b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
444b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
445b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
44642c57489SHong Zhang   merge->buf_ri    = buf_ri;
44742c57489SHong Zhang   merge->buf_rj    = buf_rj;
448de0260b3SHong Zhang   merge->owners_co = owners_co;
44977584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
45077584fe6SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
451cc6cb787SHong Zhang 
45277584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
453a914927fSHong Zhang   Cmpi->assembled      = PETSC_FALSE;
45477584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
45577584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
45642c57489SHong Zhang 
45777584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
45877584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
459f8487c73SHong Zhang   c->ptap     = ptap;
46077584fe6SHong Zhang   ptap->api   = api;
46177584fe6SHong Zhang   ptap->apj   = apj;
462d6ab1dc0SHong Zhang   ptap->rmax  = ap_rmax;
46377584fe6SHong Zhang   *C          = Cmpi;
464414894bdSHong Zhang 
465414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
466414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
467414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
468414894bdSHong Zhang   /* set default scalable */
469414894bdSHong Zhang   ptap->scalable = PETSC_TRUE;
4702205254eSKarl Rupp 
4710298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
472414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
4731795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
474414894bdSHong Zhang   } else {
4751795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
476414894bdSHong Zhang   }
477414894bdSHong Zhang 
478f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
47924ecddacSHong Zhang   if (pti[pn] != 0) {
48057622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
48157622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
482f2b054eeSHong Zhang   } else {
48377584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
484f2b054eeSHong Zhang   }
485f2b054eeSHong Zhang #endif
48642c57489SHong Zhang   PetscFunctionReturn(0);
48742c57489SHong Zhang }
48842c57489SHong Zhang 
48942c57489SHong Zhang #undef __FUNCT__
49042c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
49142c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
49242c57489SHong Zhang {
49342c57489SHong Zhang   PetscErrorCode      ierr;
494f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->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;
498f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
4999f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
5001d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
50182412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
50282412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
503e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
504d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
505ce94432eSBarry Smith   MPI_Comm            comm;
50642c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
50782412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
50842c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
50950f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
51042c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
51142c57489SHong Zhang   MPI_Status          *status;
51282412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
51382412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
514d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
5156ce94e8fSHong Zhang   PetscBool           scalable;
51638571addSHong Zhang #if defined(PTAP_PROFILE)
517b091e509SHong Zhang   PetscLogDouble t0,t1,t2,t3,t4,et2_AP=0.0,et2_PtAP=0.0,t2_0,t2_1,t2_2;
51838571addSHong Zhang #endif
51942c57489SHong Zhang 
52042c57489SHong Zhang   PetscFunctionBegin;
521ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
52238571addSHong Zhang #if defined(PTAP_PROFILE)
5238563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
52438571addSHong Zhang #endif
52542c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
52642c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
52742c57489SHong Zhang 
528f8487c73SHong Zhang   ptap = c->ptap;
529ce94432eSBarry Smith   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
530f8487c73SHong Zhang   merge    = ptap->merge;
531414894bdSHong Zhang   apa      = ptap->apa;
5326ce94e8fSHong Zhang   scalable = ptap->scalable;
5336c6e00e0SHong Zhang 
534a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
535*6b911d16SHong Zhang   /*-----------------------------------------------------*/
536f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
5379f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
538f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5396c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
540b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
541a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5426c6e00e0SHong Zhang   }
54338571addSHong Zhang #if defined(PTAP_PROFILE)
5448563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
54538571addSHong Zhang #endif
546f8487c73SHong Zhang 
547f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
548f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
54942c57489SHong Zhang   /* get data from symbolic products */
550a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
551a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
552a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
55342c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
55442c57489SHong Zhang 
555de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
5561795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
557de0260b3SHong Zhang 
558de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5597a2fc3feSBarry Smith   owners = merge->rowmap->range;
5601795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
56142c57489SHong Zhang 
562a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
563d9824c15SHong Zhang 
56438571addSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but fast */
565b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
56638571addSHong Zhang     /*-----------------------------------------------------------------------------------------------------*/
567a9d06231SHong Zhang     for (i=0; i<am; i++) {
568b091e509SHong Zhang #if defined(PTAP_PROFILE)
5698563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
570b091e509SHong Zhang #endif
571a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
572a9d06231SHong Zhang       /*------------------------------------------------------------*/
573a9d06231SHong Zhang       apJ = apj + api[i];
574a9d06231SHong Zhang 
575a9d06231SHong Zhang       /* diagonal portion of A */
576a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
577a9d06231SHong Zhang       adj = ad->j + adi[i];
578a9d06231SHong Zhang       ada = ad->a + adi[i];
579a9d06231SHong Zhang       for (j=0; j<anz; j++) {
580a9d06231SHong Zhang         row = adj[j];
581a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
582a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
583a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
584a9d06231SHong Zhang 
585a9d06231SHong Zhang         /* perform dense axpy */
586e900a757SHong Zhang         valtmp = ada[j];
587a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
588e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
589a9d06231SHong Zhang         }
590a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
591a9d06231SHong Zhang       }
592a9d06231SHong Zhang 
593a9d06231SHong Zhang       /* off-diagonal portion of A */
594a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
595a9d06231SHong Zhang       aoj = ao->j + aoi[i];
596a9d06231SHong Zhang       aoa = ao->a + aoi[i];
597a9d06231SHong Zhang       for (j=0; j<anz; j++) {
598a9d06231SHong Zhang         row = aoj[j];
599a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
600a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
601a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
602a9d06231SHong Zhang 
603a9d06231SHong Zhang         /* perform dense axpy */
604e900a757SHong Zhang         valtmp = aoa[j];
605a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
606e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
607a9d06231SHong Zhang         }
608a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
609a9d06231SHong Zhang       }
610b091e509SHong Zhang #if defined(PTAP_PROFILE)
6118563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
612b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
613b091e509SHong Zhang #endif
614a9d06231SHong Zhang 
615a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
616a9d06231SHong Zhang       /*--------------------------------------------------------------*/
617a9d06231SHong Zhang       apnz = api[i+1] - api[i];
618a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
619a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
620e900a757SHong Zhang       poJ = po->j + po->i[i];
621a9d06231SHong Zhang       pA  = po->a + po->i[i];
622a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
623e900a757SHong Zhang         row = poJ[j];
624e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
625e900a757SHong Zhang         cj  = coj + coi[row];
626e900a757SHong Zhang         ca  = coa + coi[row];
627a9d06231SHong Zhang         /* perform dense axpy */
628e900a757SHong Zhang         valtmp = pA[j];
629a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
630e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
631a9d06231SHong Zhang         }
632a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
633a9d06231SHong Zhang       }
634a9d06231SHong Zhang 
635a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
636a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
637e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
638a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
639a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
640e900a757SHong Zhang         row = pdJ[j];
641e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
642e900a757SHong Zhang         cj  = bj + bi[row];
643e900a757SHong Zhang         ca  = ba + bi[row];
644a9d06231SHong Zhang         /* perform dense axpy */
645e900a757SHong Zhang         valtmp = pA[j];
646a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
647e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
648a9d06231SHong Zhang         }
649a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
650a9d06231SHong Zhang       }
651a9d06231SHong Zhang 
652a9d06231SHong Zhang       /* zero the current row of A*P */
653a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
654b091e509SHong Zhang #if defined(PTAP_PROFILE)
6558563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
656b091e509SHong Zhang       et2_PtAP += t2_2 - t2_1;
657b091e509SHong Zhang #endif
658a9d06231SHong Zhang     }
65938571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
660b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
66138571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
662a9d06231SHong Zhang     pA=pa_loc;
66342c57489SHong Zhang     for (i=0; i<am; i++) {
664b091e509SHong Zhang #if defined(PTAP_PROFILE)
6658563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
666b091e509SHong Zhang #endif
667f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
66882412ba7SHong Zhang       apnz = api[i+1] - api[i];
66982412ba7SHong Zhang       apJ  = apj + api[i];
67042c57489SHong Zhang       /* diagonal portion of A */
67182412ba7SHong Zhang       anz = adi[i+1] - adi[i];
6721d633602SHong Zhang       adj = ad->j + adi[i];
6731d633602SHong Zhang       ada = ad->a + adi[i];
67482412ba7SHong Zhang       for (j=0; j<anz; j++) {
6751d633602SHong Zhang         row    = adj[j];
67682412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
67782412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
67882412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
679e900a757SHong Zhang         valtmp = ada[j];
680f4a743e1SHong Zhang         nextp  = 0;
68182412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
68282412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
683e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
68442c57489SHong Zhang           }
68542c57489SHong Zhang         }
686dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
68742c57489SHong Zhang       }
68842c57489SHong Zhang       /* off-diagonal portion of A */
68982412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
6901d633602SHong Zhang       aoj = ao->j + aoi[i];
6911d633602SHong Zhang       aoa = ao->a + aoi[i];
69282412ba7SHong Zhang       for (j=0; j<anz; j++) {
6931d633602SHong Zhang         row    = aoj[j];
69482412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
69582412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
69682412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
697e900a757SHong Zhang         valtmp = aoa[j];
698f4a743e1SHong Zhang         nextp  = 0;
69982412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
70082412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
701e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
70242c57489SHong Zhang           }
70342c57489SHong Zhang         }
704dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
70542c57489SHong Zhang       }
706b091e509SHong Zhang #if defined(PTAP_PROFILE)
7078563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
708b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
709b091e509SHong Zhang #endif
71042c57489SHong Zhang 
711a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
712a9d06231SHong Zhang       /*--------------------------------------------------------------*/
71382412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
714f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
71582412ba7SHong Zhang       for (j=0; j<pnz; j++) {
71642c57489SHong Zhang         nextap = 0;
717f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
71882412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
719e900a757SHong Zhang           row = *poJ;
720e900a757SHong Zhang           cj  = coj + coi[row];
721e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
72282412ba7SHong Zhang         } else {                            /* put the value into Cd */
723e900a757SHong Zhang           row = *pdJ;
724e900a757SHong Zhang           cj  = bj + bi[row];
725e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
72642c57489SHong Zhang         }
727e900a757SHong Zhang         valtmp = pA[j];
72882412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
729e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
73042c57489SHong Zhang         }
731dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
732de0260b3SHong Zhang       }
7330496f32cSHong Zhang       pA += pnz;
73442c57489SHong Zhang       /* zero the current row info for A*P */
73582412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
736b091e509SHong Zhang #if defined(PTAP_PROFILE)
7378563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
738b091e509SHong Zhang       et2_PtAP += t2_2 - t2_1;
739b091e509SHong Zhang #endif
74042c57489SHong Zhang     }
741d9824c15SHong Zhang   }
74238571addSHong Zhang #if defined(PTAP_PROFILE)
7438563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
74438571addSHong Zhang #endif
74542c57489SHong Zhang 
746a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
747a9d06231SHong Zhang   /*------------------------------------*/
74842c57489SHong Zhang   buf_ri = merge->buf_ri;
74942c57489SHong Zhang   buf_rj = merge->buf_rj;
75042c57489SHong Zhang   len_s  = merge->len_s;
751fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
75242c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
75342c57489SHong Zhang 
754dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
75542c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
75642c57489SHong Zhang     if (!len_s[proc]) continue;
757de0260b3SHong Zhang     i    = merge->owners_co[proc];
758de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
75942c57489SHong Zhang     k++;
76042c57489SHong Zhang   }
7610c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
7620c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
76342c57489SHong Zhang 
764a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
76542c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
76682412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
76738571addSHong Zhang #if defined(PTAP_PROFILE)
7688563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
76938571addSHong Zhang #endif
77042c57489SHong Zhang 
771a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
772a9d06231SHong Zhang   /*------------------------------------------------------*/
773dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
77442c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
77542c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
77642c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
77742c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
77882412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
77942c57489SHong Zhang   }
78042c57489SHong Zhang 
781de0260b3SHong Zhang   for (i=0; i<cm; i++) {
78282412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
78342c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
784de0260b3SHong Zhang     ba_i = ba + bi[i];
78582412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
78642c57489SHong Zhang     /* add received vals into ba */
78742c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
78842c57489SHong Zhang       /* i-th row */
78942c57489SHong Zhang       if (i == *nextrow[k]) {
79082412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
79182412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
79282412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
79382412ba7SHong Zhang         nextcj = 0;
79482412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
79582412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
79682412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
79742c57489SHong Zhang           }
79842c57489SHong Zhang         }
79982412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
800c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
80142c57489SHong Zhang       }
80242c57489SHong Zhang     }
80382412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
80442c57489SHong Zhang   }
80542c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80642c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80742c57489SHong Zhang 
808de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
809c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
81042c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
8110572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
81238571addSHong Zhang #if defined(PTAP_PROFILE)
8138563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
814b091e509SHong Zhang   if (rank==1) PetscPrintf(MPI_COMM_SELF,"  [%d] PtAPNum %g/P + %g/PtAP( %g + %g ) + %g/comm + %g/Cloc = %g\n\n",rank,t1-t0,t2-t1,et2_AP,et2_PtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr);
81538571addSHong Zhang #endif
81642c57489SHong Zhang   PetscFunctionReturn(0);
81742c57489SHong Zhang }
818