xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision a914927fa2c1e628be18df710bc2fdd407d2e546)
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 
1209573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
13cc6cb787SHong Zhang #undef __FUNCT__
14f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
15f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
16cc6cb787SHong Zhang {
17cc6cb787SHong Zhang   PetscErrorCode       ierr;
18f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
19f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
20cc6cb787SHong Zhang 
21cc6cb787SHong Zhang   PetscFunctionBegin;
22f8487c73SHong Zhang   if (ptap){
23c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
24b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
25f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
26a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
27a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
28c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
29c5af039cSHong Zhang     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
30c5af039cSHong Zhang     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
31414894bdSHong Zhang     if (ptap->apa){ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
32c8b0795cSMark F. Adams     if (merge) {
33cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
34cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
35cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
36cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
37cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
38c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
39cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
40c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
41cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4205b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4305b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4405b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
456bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
46dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
47f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
48bf0cc555SLisandro Dalcin     }
49f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
50c8b0795cSMark F. Adams   }
51cc6cb787SHong Zhang   PetscFunctionReturn(0);
52cc6cb787SHong Zhang }
53cc6cb787SHong Zhang 
54cc6cb787SHong Zhang #undef __FUNCT__
55cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
56f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
57f0c0a3a6SBarry Smith {
58cc6cb787SHong Zhang   PetscErrorCode       ierr;
5911a6856fSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
6011a6856fSHong Zhang   Mat_PtAPMPI          *ptap = a->ptap;
6111a6856fSHong Zhang   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
62cc6cb787SHong Zhang 
63cc6cb787SHong Zhang   PetscFunctionBegin;
64dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
6511a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
6611a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
67cc6cb787SHong Zhang   PetscFunctionReturn(0);
68cc6cb787SHong Zhang }
69cc6cb787SHong Zhang 
7042c57489SHong Zhang #undef __FUNCT__
71cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
72cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
7342c57489SHong Zhang {
7442c57489SHong Zhang   PetscErrorCode ierr;
7542c57489SHong Zhang 
7642c57489SHong Zhang   PetscFunctionBegin;
77cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
78cf3ca8ceSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
79cf3ca8ceSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
80cf3ca8ceSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
817a7894deSKris Buschelman   }
82cf3ca8ceSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
83cf3ca8ceSHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
84cf3ca8ceSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
8542c57489SHong Zhang   PetscFunctionReturn(0);
8642c57489SHong Zhang }
8742c57489SHong Zhang 
8842c57489SHong Zhang #undef __FUNCT__
8942c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9042c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9142c57489SHong Zhang {
9242c57489SHong Zhang   PetscErrorCode       ierr;
9377584fe6SHong Zhang   Mat                  Cmpi;
94f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
95a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
96f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
9742c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
9842c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
99ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10077584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
101*a914927fSHong Zhang   PetscInt             *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
102d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
10342c57489SHong Zhang   PetscBT              lnkbt;
1047adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
105*a914927fSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
10642c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
10742c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
108ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
10942c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
110ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
111ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11242c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
11377584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
114d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
115*a914927fSHong Zhang   PetscInt             rmax;
11642c57489SHong Zhang 
11742c57489SHong Zhang   PetscFunctionBegin;
118c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
119c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
120c5af039cSHong 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);
121c5af039cSHong Zhang   }
122c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
123c5af039cSHong 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);
124c5af039cSHong Zhang   }
125c5af039cSHong Zhang 
12683445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
12783445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12883445d95SHong Zhang 
129f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
130f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
131f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
1326c6e00e0SHong Zhang 
1336c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
134b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
135fe615159SHong Zhang 
1366c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
137a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1386c6e00e0SHong Zhang 
139a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
140a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1416c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1426c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
14342c57489SHong Zhang 
144fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
145fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
14677584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
14782412ba7SHong Zhang   api[0] = 0;
14883445d95SHong Zhang 
14983445d95SHong Zhang   /* create and initialize a linked list */
150*a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
15183445d95SHong Zhang 
152*a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
153d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
154f4a743e1SHong Zhang   current_space = free_space;
155f4a743e1SHong Zhang 
156f4a743e1SHong Zhang   for (i=0; i<am; i++) {
157f4a743e1SHong Zhang     /* diagonal portion of A */
158ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
15977584fe6SHong Zhang     aj  = ad->j + adi[i];
160ded4f561SHong Zhang     for (j=0; j<nzi; j++){
16177584fe6SHong Zhang       row  = aj[j];
16282412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
163ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
16483445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
165*a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
166f4a743e1SHong Zhang     }
167f4a743e1SHong Zhang     /* off-diagonal portion of A */
168ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
16977584fe6SHong Zhang     aj  = ao->j + aoi[i];
170ded4f561SHong Zhang     for (j=0; j<nzi; j++){
17177584fe6SHong Zhang       row = aj[j];
17282412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
173ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
174*a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
175f4a743e1SHong Zhang     }
176*a914927fSHong Zhang     apnz = lnk[0];
17782412ba7SHong Zhang     api[i+1] = api[i] + apnz;
17877584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
179f4a743e1SHong Zhang 
18083445d95SHong Zhang     /* if free space is not available, double the total space in the list */
18182412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1822ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
183f2b054eeSHong Zhang       nspacedouble++;
184f4a743e1SHong Zhang     }
185f4a743e1SHong Zhang 
186f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
187*a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
18882412ba7SHong Zhang     current_space->array           += apnz;
18982412ba7SHong Zhang     current_space->local_used      += apnz;
19082412ba7SHong Zhang     current_space->local_remaining -= apnz;
191f4a743e1SHong Zhang   }
192*a914927fSHong Zhang 
19382412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
194f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
19577584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
19677584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
197118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
198d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
199f4a743e1SHong Zhang 
200ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
201ded4f561SHong Zhang   /*----------------------------------------------------*/
202de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
20342c57489SHong Zhang 
204ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
205d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
20683445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
207de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
208de0260b3SHong Zhang   coi[0] = 0;
20942c57489SHong Zhang 
210d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
211d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
212a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
21342c57489SHong Zhang   current_space = free_space;
21442c57489SHong Zhang 
215de0260b3SHong Zhang   for (i=0; i<pon; i++) {
21682412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
21777584fe6SHong Zhang     ptJ = potj + poti[i];
21877584fe6SHong Zhang     for (j=0; j<pnz; j++){
21977584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
22082412ba7SHong Zhang       apnz = api[row+1] - api[row];
221ded4f561SHong Zhang       Jptr = apj + api[row];
22283445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
223*a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
22442c57489SHong Zhang     }
225*a914927fSHong Zhang     nnz = lnk[0];
22642c57489SHong Zhang 
22783445d95SHong Zhang     /* If free space is not available, double the total space in the list */
228ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2294238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
230d16ca5f0SHong Zhang       nspacedouble++;
23142c57489SHong Zhang     }
23242c57489SHong Zhang 
23342c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
234*a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
235ded4f561SHong Zhang     current_space->array           += nnz;
236ded4f561SHong Zhang     current_space->local_used      += nnz;
237ded4f561SHong Zhang     current_space->local_remaining -= nnz;
238ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
23942c57489SHong Zhang   }
240de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
241a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
242118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
243d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
244de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
24542c57489SHong Zhang 
246ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
247ded4f561SHong Zhang   /*----------------------------------------------*/
24842c57489SHong Zhang   /* determine row ownership */
24983445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
25026283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2517a2fc3feSBarry Smith   merge->rowmap->n  = pn;
2527a2fc3feSBarry Smith   merge->rowmap->bs = 1;
25326283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2547a2fc3feSBarry Smith   owners = merge->rowmap->range;
25542c57489SHong Zhang 
25642c57489SHong Zhang   /* determine the number of messages to send, their lengths */
25742c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
25883445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
25942c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
26042c57489SHong Zhang   len_s = merge->len_s;
26142c57489SHong Zhang   merge->nsend = 0;
26283445d95SHong Zhang 
263de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
264de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
265de0260b3SHong Zhang 
26683445d95SHong Zhang   proc = 0;
267de0260b3SHong Zhang   for (i=0; i<pon; i++){
268de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
269de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
270de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
27183445d95SHong Zhang   }
272de0260b3SHong Zhang 
273de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
274de0260b3SHong Zhang   owners_co[0] = 0;
27542c57489SHong Zhang   for (proc=0; proc<size; proc++){
276de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
27783445d95SHong Zhang     if (len_si[proc]){
27842c57489SHong Zhang       merge->nsend++;
27983445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
28042c57489SHong Zhang       len += len_si[proc];
28142c57489SHong Zhang     }
28242c57489SHong Zhang   }
28342c57489SHong Zhang 
284ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
28542c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
28642c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
28742c57489SHong Zhang 
288ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
289529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
290529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
291ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
29242c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
29342c57489SHong Zhang     if (!len_s[proc]) continue;
294de0260b3SHong Zhang     i = owners_co[proc];
295529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
29642c57489SHong Zhang     k++;
29742c57489SHong Zhang   }
29842c57489SHong Zhang 
299ded4f561SHong Zhang   /* receives and sends of coj are complete */
300ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
30177584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
302c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
303ded4f561SHong Zhang   }
304ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3050c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
30682412ba7SHong Zhang 
307ded4f561SHong Zhang   /* send and recv coi */
308ded4f561SHong Zhang   /*-------------------*/
309529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
310529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
31142c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
31242c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
31342c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
31442c57489SHong Zhang     if (!len_s[proc]) continue;
31542c57489SHong Zhang     /* form outgoing message for i-structure:
31642c57489SHong Zhang          buf_si[0]:                 nrows to be sent
31742c57489SHong Zhang                [1:nrows]:           row index (global)
31842c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
31942c57489SHong Zhang     */
32042c57489SHong Zhang     /*-------------------------------------------*/
32142c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
32242c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
32342c57489SHong Zhang     buf_si[0]   = nrows;
32442c57489SHong Zhang     buf_si_i[0] = 0;
32542c57489SHong Zhang     nrows = 0;
326de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
327ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
328ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
329de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
33042c57489SHong Zhang       nrows++;
33142c57489SHong Zhang     }
332529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
33342c57489SHong Zhang     k++;
33442c57489SHong Zhang     buf_si += len_si[proc];
33542c57489SHong Zhang   }
336ded4f561SHong Zhang   i = merge->nrecv;
337ded4f561SHong Zhang   while (i--) {
338c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
339ded4f561SHong Zhang   }
340ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3410c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
342*a914927fSHong Zhang 
34342c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
34442c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
345ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
346ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
34742c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
34842c57489SHong Zhang 
349de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
350de0260b3SHong Zhang   /*------------------------------------------*/
351ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
352ded4f561SHong Zhang 
35342c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
35442c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
35542c57489SHong Zhang   bi[0] = 0;
35642c57489SHong Zhang 
357d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
358d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
359a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
36042c57489SHong Zhang   current_space = free_space;
36142c57489SHong Zhang 
362687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
36342c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
36442c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
36542c57489SHong Zhang     nrows       = *buf_ri_k[k];
36642c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
36742c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
36842c57489SHong Zhang   }
36942c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
37008cb4532SHong Zhang   rmax = 0;
37142c57489SHong Zhang   for (i=0; i<pn; i++) {
372ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
373ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
37477584fe6SHong Zhang     ptJ = pdtj + pdti[i];
37577584fe6SHong Zhang     for (j=0; j<pnz; j++){
37677584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
377ded4f561SHong Zhang       apnz = api[row+1] - api[row];
378ded4f561SHong Zhang       Jptr = apj + api[row];
379ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
380*a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
381ded4f561SHong Zhang     }
382*a914927fSHong Zhang 
38342c57489SHong Zhang     /* add received col data into lnk */
38442c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
38542c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
386ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
387ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
388*a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
38942c57489SHong Zhang         nextrow[k]++; nextci[k]++;
39042c57489SHong Zhang       }
39142c57489SHong Zhang     }
392*a914927fSHong Zhang     nnz = lnk[0];
39342c57489SHong Zhang 
39442c57489SHong Zhang     /* if free space is not available, make more free space */
395ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
3964238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
397d16ca5f0SHong Zhang       nspacedouble++;
39842c57489SHong Zhang     }
39942c57489SHong Zhang     /* copy data into free space, then initialize lnk */
400*a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
401ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
402ded4f561SHong Zhang     current_space->array           += nnz;
403ded4f561SHong Zhang     current_space->local_used      += nnz;
404ded4f561SHong Zhang     current_space->local_remaining -= nnz;
405ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
40608cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
40742c57489SHong Zhang   }
408ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4090572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
41042c57489SHong Zhang 
41142c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
412a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
413118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]+1);
414d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
41542c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
41642c57489SHong Zhang 
417*a914927fSHong Zhang   /* create symbolic parallel matrix Cmpi */
418*a914927fSHong Zhang   /*--------------------------------------*/
41977584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
42077584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
421a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
42277584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
42377584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
42442c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
425a2f3521dSMark F. Adams 
42642c57489SHong Zhang   merge->bi            = bi;
42742c57489SHong Zhang   merge->bj            = bj;
428de0260b3SHong Zhang   merge->coi           = coi;
429de0260b3SHong Zhang   merge->coj           = coj;
43042c57489SHong Zhang   merge->buf_ri        = buf_ri;
43142c57489SHong Zhang   merge->buf_rj        = buf_rj;
432de0260b3SHong Zhang   merge->owners_co     = owners_co;
43377584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
43477584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
435cc6cb787SHong Zhang 
43677584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
437*a914927fSHong Zhang   Cmpi->assembled      = PETSC_FALSE;
43877584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
43977584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
44042c57489SHong Zhang 
44177584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
44277584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
443f8487c73SHong Zhang   c->ptap        = ptap;
44477584fe6SHong Zhang   ptap->api      = api;
44577584fe6SHong Zhang   ptap->apj      = apj;
446f8487c73SHong Zhang   ptap->merge    = merge;
447d6ab1dc0SHong Zhang   ptap->rmax     = ap_rmax;
44877584fe6SHong Zhang   *C             = Cmpi;
449414894bdSHong Zhang 
450414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
451414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
452414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
453414894bdSHong Zhang   /* set default scalable */
454414894bdSHong Zhang   ptap->scalable = PETSC_TRUE;
455414894bdSHong Zhang   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,PETSC_NULL);CHKERRQ(ierr);
456414894bdSHong Zhang   if (!ptap->scalable){  /* Do dense axpy */
457414894bdSHong Zhang     ierr = PetscMalloc(pN*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
458414894bdSHong Zhang   } else {
459414894bdSHong Zhang     ierr = PetscMalloc((ap_rmax+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
460414894bdSHong Zhang   }
461414894bdSHong Zhang 
462f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
463f2b054eeSHong Zhang   if (bi[pn] != 0) {
46477584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
46577584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
466f2b054eeSHong Zhang   } else {
46777584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
468f2b054eeSHong Zhang   }
469f2b054eeSHong Zhang #endif
47042c57489SHong Zhang   PetscFunctionReturn(0);
47142c57489SHong Zhang }
47242c57489SHong Zhang 
47342c57489SHong Zhang #undef __FUNCT__
47442c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
47542c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
47642c57489SHong Zhang {
47742c57489SHong Zhang   PetscErrorCode       ierr;
47842c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
479f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
48042c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
481de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
48242c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
483f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
4841d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
48582412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
48682412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
487e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
488d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
4897adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
49042c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
49182412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
49242c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
49350f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
49442c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
49542c57489SHong Zhang   MPI_Status           *status;
49682412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
49782412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
498d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
4996ce94e8fSHong Zhang   PetscBool            scalable;
50042c57489SHong Zhang 
50142c57489SHong Zhang   PetscFunctionBegin;
50242c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
50342c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
50442c57489SHong Zhang 
505f8487c73SHong Zhang   ptap = c->ptap;
5061c4805f8SJed 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");
507f8487c73SHong Zhang   merge    = ptap->merge;
508414894bdSHong Zhang   apa      = ptap->apa;
5096ce94e8fSHong Zhang   scalable = ptap->scalable;
5106c6e00e0SHong Zhang 
511a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
512f9473cf7SHong Zhang   /*--------------------------------------------------*/
513f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
514f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5156c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
516b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
517a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5186c6e00e0SHong Zhang   }
519f8487c73SHong Zhang 
520f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
521f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
52242c57489SHong Zhang   /* get data from symbolic products */
523a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
524a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
525a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
52642c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
52742c57489SHong Zhang 
528de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
529de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
530de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
531de0260b3SHong Zhang 
532de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5337a2fc3feSBarry Smith   owners = merge->rowmap->range;
534de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
535de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
53642c57489SHong Zhang 
537a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
538d9824c15SHong Zhang 
539d9824c15SHong Zhang   if (!scalable){  /* Do dense axpy */
540a9d06231SHong Zhang     /*--------------------------------------------------*/
541414894bdSHong Zhang     /* apa (length of pN) stores dense row A[i,:]*P - nonscalable! */
542a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
543a9d06231SHong Zhang 
544a9d06231SHong Zhang     for (i=0; i<am; i++) {
545a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
546a9d06231SHong Zhang       /*------------------------------------------------------------*/
547a9d06231SHong Zhang       apJ = apj + api[i];
548a9d06231SHong Zhang 
549a9d06231SHong Zhang       /* diagonal portion of A */
550a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
551a9d06231SHong Zhang       adj = ad->j + adi[i];
552a9d06231SHong Zhang       ada = ad->a + adi[i];
553a9d06231SHong Zhang       for (j=0; j<anz; j++) {
554a9d06231SHong Zhang         row = adj[j];
555a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
556a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
557a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
558a9d06231SHong Zhang 
559a9d06231SHong Zhang         /* perform dense axpy */
560e900a757SHong Zhang         valtmp = ada[j];
561a9d06231SHong Zhang         for (k=0; k<pnz; k++){
562e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
563a9d06231SHong Zhang         }
564a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
565a9d06231SHong Zhang       }
566a9d06231SHong Zhang 
567a9d06231SHong Zhang       /* off-diagonal portion of A */
568a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
569a9d06231SHong Zhang       aoj = ao->j + aoi[i];
570a9d06231SHong Zhang       aoa = ao->a + aoi[i];
571a9d06231SHong Zhang       for (j=0; j<anz; j++) {
572a9d06231SHong Zhang         row = aoj[j];
573a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
574a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
575a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
576a9d06231SHong Zhang 
577a9d06231SHong Zhang         /* perform dense axpy */
578e900a757SHong Zhang         valtmp = aoa[j];
579a9d06231SHong Zhang         for (k=0; k<pnz; k++){
580e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
581a9d06231SHong Zhang         }
582a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
583a9d06231SHong Zhang       }
584a9d06231SHong Zhang 
585a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
586a9d06231SHong Zhang       /*--------------------------------------------------------------*/
587a9d06231SHong Zhang       apnz = api[i+1] - api[i];
588a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
589a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
590e900a757SHong Zhang       poJ = po->j + po->i[i];
591a9d06231SHong Zhang       pA  = po->a + po->i[i];
592a9d06231SHong Zhang       for (j=0; j<pnz; j++){
593e900a757SHong Zhang         row = poJ[j];
594e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
595e900a757SHong Zhang         cj  = coj + coi[row];
596e900a757SHong Zhang         ca  = coa + coi[row];
597a9d06231SHong Zhang         /* perform dense axpy */
598e900a757SHong Zhang         valtmp = pA[j];
599a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
600e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
601a9d06231SHong Zhang         }
602a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
603a9d06231SHong Zhang       }
604a9d06231SHong Zhang 
605a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
606a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
607e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
608a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
609a9d06231SHong Zhang       for (j=0; j<pnz; j++){
610e900a757SHong Zhang         row = pdJ[j];
611e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
612e900a757SHong Zhang         cj  = bj + bi[row];
613e900a757SHong Zhang         ca  = ba + bi[row];
614a9d06231SHong Zhang         /* perform dense axpy */
615e900a757SHong Zhang         valtmp = pA[j];
616a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
617e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
618a9d06231SHong Zhang         }
619a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
620a9d06231SHong Zhang       }
621a9d06231SHong Zhang 
622a9d06231SHong Zhang       /* zero the current row of A*P */
623a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
624a9d06231SHong Zhang     }
625d9824c15SHong Zhang   } else {/* Perform sparse axpy */
626a9d06231SHong Zhang     /*----------------------------------------------------*/
627414894bdSHong Zhang     /* apa (length ap_rmax) stores sparse row A[i,:]*P */
628d6ab1dc0SHong Zhang     ierr = PetscMemzero(apa,ptap->rmax*sizeof(MatScalar));CHKERRQ(ierr);
62942c57489SHong Zhang 
630a9d06231SHong Zhang     pA=pa_loc;
63142c57489SHong Zhang     for (i=0; i<am; i++) {
632f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
63382412ba7SHong Zhang       apnz = api[i+1] - api[i];
63482412ba7SHong Zhang       apJ  = apj + api[i];
63542c57489SHong Zhang       /* diagonal portion of A */
63682412ba7SHong Zhang       anz = adi[i+1] - adi[i];
6371d633602SHong Zhang       adj = ad->j + adi[i];
6381d633602SHong Zhang       ada = ad->a + adi[i];
63982412ba7SHong Zhang       for (j=0; j<anz; j++) {
6401d633602SHong Zhang         row = adj[j];
64182412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
64282412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
64382412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
644e900a757SHong Zhang         valtmp = ada[j];
645f4a743e1SHong Zhang         nextp  = 0;
64682412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
64782412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
648e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
64942c57489SHong Zhang           }
65042c57489SHong Zhang         }
651dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
65242c57489SHong Zhang       }
65342c57489SHong Zhang       /* off-diagonal portion of A */
65482412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
6551d633602SHong Zhang       aoj = ao->j + aoi[i];
6561d633602SHong Zhang       aoa = ao->a + aoi[i];
65782412ba7SHong Zhang       for (j=0; j<anz; j++) {
6581d633602SHong Zhang         row = aoj[j];
65982412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
66082412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
66182412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
662e900a757SHong Zhang         valtmp = aoa[j];
663f4a743e1SHong Zhang         nextp  = 0;
66482412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
66582412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
666e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
66742c57489SHong Zhang           }
66842c57489SHong Zhang         }
669dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
67042c57489SHong Zhang       }
67142c57489SHong Zhang 
672a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
673a9d06231SHong Zhang       /*--------------------------------------------------------------*/
67482412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
675f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
67682412ba7SHong Zhang       for (j=0; j<pnz; j++) {
67742c57489SHong Zhang         nextap = 0;
678f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
67982412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
680e900a757SHong Zhang           row = *poJ;
681e900a757SHong Zhang           cj  = coj + coi[row];
682e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
68382412ba7SHong Zhang         } else {                            /* put the value into Cd */
684e900a757SHong Zhang           row  = *pdJ;
685e900a757SHong Zhang           cj   = bj + bi[row];
686e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
68742c57489SHong Zhang         }
688e900a757SHong Zhang         valtmp = pA[j];
68982412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
690e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
69142c57489SHong Zhang         }
692dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
693de0260b3SHong Zhang       }
6940496f32cSHong Zhang       pA += pnz;
69542c57489SHong Zhang       /* zero the current row info for A*P */
69682412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
69742c57489SHong Zhang     }
698d9824c15SHong Zhang   }
69942c57489SHong Zhang 
700a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
701a9d06231SHong Zhang   /*------------------------------------*/
70242c57489SHong Zhang   buf_ri = merge->buf_ri;
70342c57489SHong Zhang   buf_rj = merge->buf_rj;
70442c57489SHong Zhang   len_s  = merge->len_s;
705fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
70642c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
70742c57489SHong Zhang 
708a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
70942c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
71042c57489SHong Zhang     if (!len_s[proc]) continue;
711de0260b3SHong Zhang     i = merge->owners_co[proc];
712de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
71342c57489SHong Zhang     k++;
71442c57489SHong Zhang   }
7150c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
7160c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
71742c57489SHong Zhang 
718a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
71942c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
72082412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
72142c57489SHong Zhang 
722a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
723a9d06231SHong Zhang   /*------------------------------------------------------*/
724687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
72542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
72642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
72742c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
72842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
72982412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
73042c57489SHong Zhang   }
73142c57489SHong Zhang 
732de0260b3SHong Zhang   for (i=0; i<cm; i++) {
73382412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
73442c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
735de0260b3SHong Zhang     ba_i = ba + bi[i];
73682412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
73742c57489SHong Zhang     /* add received vals into ba */
73842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
73942c57489SHong Zhang       /* i-th row */
74042c57489SHong Zhang       if (i == *nextrow[k]) {
74182412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
74282412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
74382412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
74482412ba7SHong Zhang         nextcj = 0;
74582412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
74682412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
74782412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
74842c57489SHong Zhang           }
74942c57489SHong Zhang         }
75082412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
751c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
75242c57489SHong Zhang       }
75342c57489SHong Zhang     }
75482412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
75542c57489SHong Zhang   }
75642c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
75742c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
75842c57489SHong Zhang 
759de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
760c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
76142c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
7620572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
763d9824c15SHong Zhang #if defined(PETSC_USE_INFO)
764d9824c15SHong Zhang   if (scalable){
765d9824c15SHong Zhang     ierr = PetscInfo(C,"Use scalable sparse axpy\n");CHKERRQ(ierr);
766d9824c15SHong Zhang   } else {
767d9824c15SHong Zhang     ierr = PetscInfo(C,"Use non-scalable dense axpy\n");CHKERRQ(ierr);
768d9824c15SHong Zhang   }
769d9824c15SHong Zhang #endif
77042c57489SHong Zhang   PetscFunctionReturn(0);
77142c57489SHong Zhang }
772