xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 6ce94e8fe7f45068660bdbda61d05c45a2a71fa4)
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 
12fe615159SHong Zhang /*
13fe615159SHong Zhang #define DEBUG_MATPTAP
14fe615159SHong Zhang  */
15fe615159SHong Zhang 
1609573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
17cc6cb787SHong Zhang #undef __FUNCT__
18f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
19f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
20cc6cb787SHong Zhang {
21cc6cb787SHong Zhang   PetscErrorCode       ierr;
22f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
23f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
24cc6cb787SHong Zhang 
25cc6cb787SHong Zhang   PetscFunctionBegin;
26f8487c73SHong Zhang   if (ptap){
27c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
28b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
29f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
30a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
31a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
32c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
33c5af039cSHong Zhang     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
34c5af039cSHong Zhang     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
35414894bdSHong Zhang     if (ptap->apa){ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
36c8b0795cSMark F. Adams     if (merge) {
37cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
38cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
39cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
40cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
41cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
42c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
43cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
44c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
45cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4605b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4705b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4805b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
496bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
50dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
51f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
52bf0cc555SLisandro Dalcin     }
53f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
54c8b0795cSMark F. Adams   }
55c8b0795cSMark F. Adams 
56cc6cb787SHong Zhang   PetscFunctionReturn(0);
57cc6cb787SHong Zhang }
58cc6cb787SHong Zhang 
59cc6cb787SHong Zhang #undef __FUNCT__
60cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
61f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
62f0c0a3a6SBarry Smith {
63cc6cb787SHong Zhang   PetscErrorCode       ierr;
6411a6856fSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
6511a6856fSHong Zhang   Mat_PtAPMPI          *ptap = a->ptap;
6611a6856fSHong Zhang   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
67cc6cb787SHong Zhang 
68cc6cb787SHong Zhang   PetscFunctionBegin;
69dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
7011a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
7111a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
72cc6cb787SHong Zhang   PetscFunctionReturn(0);
73cc6cb787SHong Zhang }
74cc6cb787SHong Zhang 
7542c57489SHong Zhang #undef __FUNCT__
76cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
77cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
7842c57489SHong Zhang {
7942c57489SHong Zhang   PetscErrorCode ierr;
8042c57489SHong Zhang 
8142c57489SHong Zhang   PetscFunctionBegin;
82cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
83cf3ca8ceSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
84cf3ca8ceSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
85cf3ca8ceSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
867a7894deSKris Buschelman   }
87cf3ca8ceSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
88cf3ca8ceSHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
89cf3ca8ceSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
9042c57489SHong Zhang   PetscFunctionReturn(0);
9142c57489SHong Zhang }
9242c57489SHong Zhang 
9342c57489SHong Zhang #undef __FUNCT__
9442c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9542c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9642c57489SHong Zhang {
9742c57489SHong Zhang   PetscErrorCode       ierr;
9877584fe6SHong Zhang   Mat                  Cmpi;
99f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
100a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
101f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10242c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10342c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
104ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10577584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
10613f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
107d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
10842c57489SHong Zhang   PetscBT              lnkbt;
1097adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
110529bc5dcSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
11142c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11242c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
113ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11442c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
115ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
116ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11742c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
11877584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
119d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
12008cb4532SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax;
12108cb4532SHong Zhang   PetscScalar          *vals;
12242c57489SHong Zhang 
12342c57489SHong Zhang   PetscFunctionBegin;
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 */
136f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
137f8487c73SHong Zhang   ptap->reuse    = MAT_INITIAL_MATRIX;
1386c6e00e0SHong Zhang 
13977584fe6SHong Zhang 
1406c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
141b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
142fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
143fe615159SHong Zhang    if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
144fe615159SHong Zhang #endif
145fe615159SHong Zhang 
1466c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
147a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
148fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
149fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
150fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done...\n",rank);
151fe615159SHong Zhang #endif
1526c6e00e0SHong Zhang 
153a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
154a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1556c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1566c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
15742c57489SHong Zhang 
158fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
159fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
16077584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
16182412ba7SHong Zhang   api[0] = 0;
16283445d95SHong Zhang 
16383445d95SHong Zhang   /* create and initialize a linked list */
16483445d95SHong Zhang   nlnk = pN+1;
165fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
166fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
167fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate(), pN %D, P_rmax %D %D; Annz %D\n",rank,pN,p_loc->rmax,p_oth->rmax,adi[am]+aoi[am]);
168fe615159SHong Zhang #endif
16983445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
170fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
171fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
172fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate() is done\n",rank);
173fe615159SHong Zhang #endif
17483445d95SHong Zhang 
175d16ca5f0SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
176d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
177f4a743e1SHong Zhang   current_space = free_space;
178f4a743e1SHong Zhang 
179f4a743e1SHong Zhang   for (i=0; i<am; i++) {
18082412ba7SHong Zhang     apnz = 0;
181f4a743e1SHong Zhang     /* diagonal portion of A */
182ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
18377584fe6SHong Zhang     aj  = ad->j + adi[i];
184ded4f561SHong Zhang     for (j=0; j<nzi; j++){
18577584fe6SHong Zhang       row  = aj[j];
18682412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
187ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
18883445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
189dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
19082412ba7SHong Zhang       apnz += nlnk;
191f4a743e1SHong Zhang     }
192f4a743e1SHong Zhang     /* off-diagonal portion of A */
193ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
19477584fe6SHong Zhang     aj  = ao->j + aoi[i];
195ded4f561SHong Zhang     for (j=0; j<nzi; j++){
19677584fe6SHong Zhang       row = aj[j];
19782412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
198ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
199dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
20082412ba7SHong Zhang       apnz += nlnk;
201f4a743e1SHong Zhang     }
202f4a743e1SHong Zhang 
20382412ba7SHong Zhang     api[i+1] = api[i] + apnz;
20477584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
205f4a743e1SHong Zhang 
20683445d95SHong Zhang     /* if free space is not available, double the total space in the list */
20782412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
2082ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
209f2b054eeSHong Zhang       nspacedouble++;
210f4a743e1SHong Zhang     }
211f4a743e1SHong Zhang 
212f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
21382412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
21482412ba7SHong Zhang     current_space->array           += apnz;
21582412ba7SHong Zhang     current_space->local_used      += apnz;
21682412ba7SHong Zhang     current_space->local_remaining -= apnz;
217f4a743e1SHong Zhang   }
21882412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
219f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
22077584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
22177584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
222118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
223d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
224f4a743e1SHong Zhang 
225ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
226ded4f561SHong Zhang   /*----------------------------------------------------*/
227de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
22842c57489SHong Zhang 
229ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
230d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
23183445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
232de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
233de0260b3SHong Zhang   coi[0] = 0;
23442c57489SHong Zhang 
235d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
236d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
237a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
23842c57489SHong Zhang   current_space = free_space;
23942c57489SHong Zhang 
240de0260b3SHong Zhang   for (i=0; i<pon; i++) {
241ded4f561SHong Zhang     nnz = 0;
24282412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
24377584fe6SHong Zhang     ptJ = potj + poti[i];
24477584fe6SHong Zhang     for (j=0; j<pnz; j++){
24577584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
24682412ba7SHong Zhang       apnz = api[row+1] - api[row];
247ded4f561SHong Zhang       Jptr = apj + api[row];
24883445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
249dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
250ded4f561SHong Zhang       nnz += nlnk;
25142c57489SHong Zhang     }
25242c57489SHong Zhang 
25383445d95SHong Zhang     /* If free space is not available, double the total space in the list */
254ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2554238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
256d16ca5f0SHong Zhang       nspacedouble++;
25742c57489SHong Zhang     }
25842c57489SHong Zhang 
25942c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
260ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
261ded4f561SHong Zhang     current_space->array           += nnz;
262ded4f561SHong Zhang     current_space->local_used      += nnz;
263ded4f561SHong Zhang     current_space->local_remaining -= nnz;
264ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
26542c57489SHong Zhang   }
266de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
267a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
268118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
269d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
270de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
27142c57489SHong Zhang 
272fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
273fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
274fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] local Co is done\n",rank);
275fe615159SHong Zhang #endif
276fe615159SHong Zhang 
277ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
278ded4f561SHong Zhang   /*----------------------------------------------*/
27942c57489SHong Zhang   /* determine row ownership */
28083445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
28126283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2827a2fc3feSBarry Smith   merge->rowmap->n = pn;
2837a2fc3feSBarry Smith   merge->rowmap->bs = 1;
28426283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2857a2fc3feSBarry Smith   owners = merge->rowmap->range;
28642c57489SHong Zhang 
28742c57489SHong Zhang   /* determine the number of messages to send, their lengths */
28842c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
28983445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
29042c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
29142c57489SHong Zhang   len_s = merge->len_s;
29242c57489SHong Zhang   merge->nsend = 0;
29383445d95SHong Zhang 
294de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
295de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
296de0260b3SHong Zhang 
29783445d95SHong Zhang   proc = 0;
298de0260b3SHong Zhang   for (i=0; i<pon; i++){
299de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
300de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
301de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
30283445d95SHong Zhang   }
303de0260b3SHong Zhang 
304de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
305de0260b3SHong Zhang   owners_co[0] = 0;
30642c57489SHong Zhang   for (proc=0; proc<size; proc++){
307de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
30883445d95SHong Zhang     if (len_si[proc]){
30942c57489SHong Zhang       merge->nsend++;
31083445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
31142c57489SHong Zhang       len += len_si[proc];
31242c57489SHong Zhang     }
31342c57489SHong Zhang   }
31442c57489SHong Zhang 
315ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
31642c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
31742c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
31842c57489SHong Zhang 
319ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
320529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
321529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
322ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
32342c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
32442c57489SHong Zhang     if (!len_s[proc]) continue;
325de0260b3SHong Zhang     i = owners_co[proc];
326529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
32742c57489SHong Zhang     k++;
32842c57489SHong Zhang   }
32942c57489SHong Zhang 
330ded4f561SHong Zhang   /* receives and sends of coj are complete */
331ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
33277584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
333c280ed6aSJed Brown     PetscMPIInt icompleted;
334c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
335ded4f561SHong Zhang   }
336ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3370c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
33882412ba7SHong Zhang 
339ded4f561SHong Zhang   /* send and recv coi */
340ded4f561SHong Zhang   /*-------------------*/
341529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
342529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
34342c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
34442c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
34542c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
34642c57489SHong Zhang     if (!len_s[proc]) continue;
34742c57489SHong Zhang     /* form outgoing message for i-structure:
34842c57489SHong Zhang          buf_si[0]:                 nrows to be sent
34942c57489SHong Zhang                [1:nrows]:           row index (global)
35042c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
35142c57489SHong Zhang     */
35242c57489SHong Zhang     /*-------------------------------------------*/
35342c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
35442c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
35542c57489SHong Zhang     buf_si[0]   = nrows;
35642c57489SHong Zhang     buf_si_i[0] = 0;
35742c57489SHong Zhang     nrows = 0;
358de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
359ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
360ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
361de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
36242c57489SHong Zhang       nrows++;
36342c57489SHong Zhang     }
364529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
36542c57489SHong Zhang     k++;
36642c57489SHong Zhang     buf_si += len_si[proc];
36742c57489SHong Zhang   }
368ded4f561SHong Zhang   i = merge->nrecv;
369ded4f561SHong Zhang   while (i--) {
370c280ed6aSJed Brown     PetscMPIInt icompleted;
371c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
372ded4f561SHong Zhang   }
373ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3740c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
375f2b054eeSHong Zhang   /*
376ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
37742c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
378ae15b995SBarry Smith     ierr = PetscInfo3(A,"recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr);
37942c57489SHong Zhang   }
380f2b054eeSHong Zhang   */
38142c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
38242c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
383ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
384ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
38542c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
38642c57489SHong Zhang 
387fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
388fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
389fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Co is done\n",rank);
390fe615159SHong Zhang #endif
391fe615159SHong Zhang 
392de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
393de0260b3SHong Zhang   /*------------------------------------------*/
394ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
395ded4f561SHong Zhang 
39642c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
39742c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
39842c57489SHong Zhang   bi[0] = 0;
39942c57489SHong Zhang 
400d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
401d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
402a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
40342c57489SHong Zhang   current_space = free_space;
40442c57489SHong Zhang 
405687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
40642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
40742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
40842c57489SHong Zhang     nrows       = *buf_ri_k[k];
40942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
41042c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
41142c57489SHong Zhang   }
41242c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
41308cb4532SHong Zhang   rmax = 0;
41442c57489SHong Zhang   for (i=0; i<pn; i++) {
415ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
416ded4f561SHong Zhang     nnz = 0;
417ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
41877584fe6SHong Zhang     ptJ = pdtj + pdti[i];
41977584fe6SHong Zhang     for (j=0; j<pnz; j++){
42077584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
421ded4f561SHong Zhang       apnz = api[row+1] - api[row];
422ded4f561SHong Zhang       Jptr = apj + api[row];
423ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
424dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
425ded4f561SHong Zhang       nnz += nlnk;
426ded4f561SHong Zhang     }
42742c57489SHong Zhang     /* add received col data into lnk */
42842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
42942c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
430ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
431ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
432dadf0e6bSHong Zhang         ierr = PetscLLAddSorted(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
433ded4f561SHong Zhang         nnz += nlnk;
43442c57489SHong Zhang         nextrow[k]++; nextci[k]++;
43542c57489SHong Zhang       }
43642c57489SHong Zhang     }
43742c57489SHong Zhang 
43842c57489SHong Zhang     /* if free space is not available, make more free space */
439ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4404238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
441d16ca5f0SHong Zhang       nspacedouble++;
44242c57489SHong Zhang     }
44342c57489SHong Zhang     /* copy data into free space, then initialize lnk */
444ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
445ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
446ded4f561SHong Zhang     current_space->array           += nnz;
447ded4f561SHong Zhang     current_space->local_used      += nnz;
448ded4f561SHong Zhang     current_space->local_remaining -= nnz;
449ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
45008cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
45142c57489SHong Zhang   }
452ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4530572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
45442c57489SHong Zhang 
45542c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
456a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
457118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]+1);
458d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
45942c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
46042c57489SHong Zhang 
46108cb4532SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part - check runex55_SA ?  */
46208cb4532SHong Zhang   /*------------------------------------------------------------------------------------------------------*/
46308cb4532SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
46408cb4532SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
46508cb4532SHong Zhang 
46677584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
46777584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
468a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
46977584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
47077584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
47142c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
472a2f3521dSMark F. Adams 
47308cb4532SHong Zhang   for (i=0; i<pn; i++){
47408cb4532SHong Zhang     row = i + rstart;
47508cb4532SHong Zhang     nnz = bi[i+1] - bi[i];
47608cb4532SHong Zhang     Jptr = bj + bi[i];
47708cb4532SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
47808cb4532SHong Zhang   }
47908cb4532SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48008cb4532SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48108cb4532SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
48242c57489SHong Zhang 
48342c57489SHong Zhang   merge->bi            = bi;
48442c57489SHong Zhang   merge->bj            = bj;
485de0260b3SHong Zhang   merge->coi           = coi;
486de0260b3SHong Zhang   merge->coj           = coj;
48742c57489SHong Zhang   merge->buf_ri        = buf_ri;
48842c57489SHong Zhang   merge->buf_rj        = buf_rj;
489de0260b3SHong Zhang   merge->owners_co     = owners_co;
49077584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
49177584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
492cc6cb787SHong Zhang 
49377584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
49408cb4532SHong Zhang   /* Cmpi->assembled      = PETSC_FALSE; */
49577584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
49677584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
49742c57489SHong Zhang 
49877584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
49977584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
500f8487c73SHong Zhang   c->ptap        = ptap;
50177584fe6SHong Zhang   ptap->api      = api;
50277584fe6SHong Zhang   ptap->apj      = apj;
503f8487c73SHong Zhang   ptap->merge    = merge;
504d6ab1dc0SHong Zhang   ptap->rmax     = ap_rmax;
505f8487c73SHong Zhang 
50677584fe6SHong Zhang   *C = Cmpi;
507414894bdSHong Zhang 
508414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
509414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
510414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
511414894bdSHong Zhang   /* set default scalable */
512414894bdSHong Zhang   ptap->scalable = PETSC_TRUE;
513414894bdSHong Zhang   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,PETSC_NULL);CHKERRQ(ierr);
514414894bdSHong Zhang   if (!ptap->scalable){  /* Do dense axpy */
515414894bdSHong Zhang     ierr = PetscMalloc(pN*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
516414894bdSHong Zhang   } else {
517414894bdSHong Zhang     ierr = PetscMalloc((ap_rmax+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
518414894bdSHong Zhang   }
519414894bdSHong Zhang 
520f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
521f2b054eeSHong Zhang   if (bi[pn] != 0) {
52277584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
52377584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
524f2b054eeSHong Zhang   } else {
52577584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
526f2b054eeSHong Zhang   }
527f2b054eeSHong Zhang #endif
52842c57489SHong Zhang   PetscFunctionReturn(0);
52942c57489SHong Zhang }
53042c57489SHong Zhang 
53142c57489SHong Zhang #undef __FUNCT__
53242c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
53342c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
53442c57489SHong Zhang {
53542c57489SHong Zhang   PetscErrorCode       ierr;
53642c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
537f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
53842c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
539de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
54042c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
541f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
5421d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
54382412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
54482412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
545e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
546d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5477adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
54842c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
54982412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
55042c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
55150f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
55242c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
55342c57489SHong Zhang   MPI_Status           *status;
55482412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
55582412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
556d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
557*6ce94e8fSHong Zhang   PetscBool            scalable;
55842c57489SHong Zhang 
55942c57489SHong Zhang   PetscFunctionBegin;
56042c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
56142c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
56242c57489SHong Zhang 
563f8487c73SHong Zhang   ptap = c->ptap;
5641c4805f8SJed 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");
565f8487c73SHong Zhang   merge    = ptap->merge;
566414894bdSHong Zhang   apa      = ptap->apa;
567*6ce94e8fSHong Zhang   scalable = ptap->scalable;
5686c6e00e0SHong Zhang 
569a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
570f9473cf7SHong Zhang   /*--------------------------------------------------*/
571f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
572f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5736c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
574b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
575a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5766c6e00e0SHong Zhang   }
577f8487c73SHong Zhang 
578f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
579f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
58042c57489SHong Zhang   /* get data from symbolic products */
581a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
582a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
583a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
58442c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
58542c57489SHong Zhang 
586de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
587de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
588de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
589de0260b3SHong Zhang 
590de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5917a2fc3feSBarry Smith   owners = merge->rowmap->range;
592de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
593de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
59442c57489SHong Zhang 
595a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
596d9824c15SHong Zhang 
597d9824c15SHong Zhang   if (!scalable){  /* Do dense axpy */
598a9d06231SHong Zhang     /*--------------------------------------------------*/
599414894bdSHong Zhang     /* apa (length of pN) stores dense row A[i,:]*P - nonscalable! */
600a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
601a9d06231SHong Zhang 
602a9d06231SHong Zhang     for (i=0; i<am; i++) {
603a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
604a9d06231SHong Zhang       /*------------------------------------------------------------*/
605a9d06231SHong Zhang       apJ = apj + api[i];
606a9d06231SHong Zhang 
607a9d06231SHong Zhang       /* diagonal portion of A */
608a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
609a9d06231SHong Zhang       adj = ad->j + adi[i];
610a9d06231SHong Zhang       ada = ad->a + adi[i];
611a9d06231SHong Zhang       for (j=0; j<anz; j++) {
612a9d06231SHong Zhang         row = adj[j];
613a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
614a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
615a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
616a9d06231SHong Zhang 
617a9d06231SHong Zhang         /* perform dense axpy */
618e900a757SHong Zhang         valtmp = ada[j];
619a9d06231SHong Zhang         for (k=0; k<pnz; k++){
620e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
621a9d06231SHong Zhang         }
622a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
623a9d06231SHong Zhang       }
624a9d06231SHong Zhang 
625a9d06231SHong Zhang       /* off-diagonal portion of A */
626a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
627a9d06231SHong Zhang       aoj = ao->j + aoi[i];
628a9d06231SHong Zhang       aoa = ao->a + aoi[i];
629a9d06231SHong Zhang       for (j=0; j<anz; j++) {
630a9d06231SHong Zhang         row = aoj[j];
631a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
632a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
633a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
634a9d06231SHong Zhang 
635a9d06231SHong Zhang         /* perform dense axpy */
636e900a757SHong Zhang         valtmp = aoa[j];
637a9d06231SHong Zhang         for (k=0; k<pnz; k++){
638e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
639a9d06231SHong Zhang         }
640a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
641a9d06231SHong Zhang       }
642a9d06231SHong Zhang 
643a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
644a9d06231SHong Zhang       /*--------------------------------------------------------------*/
645a9d06231SHong Zhang       apnz = api[i+1] - api[i];
646a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
647a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
648e900a757SHong Zhang       poJ = po->j + po->i[i];
649a9d06231SHong Zhang       pA  = po->a + po->i[i];
650a9d06231SHong Zhang       for (j=0; j<pnz; j++){
651e900a757SHong Zhang         row = poJ[j];
652e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
653e900a757SHong Zhang         cj  = coj + coi[row];
654e900a757SHong Zhang         ca  = coa + coi[row];
655a9d06231SHong Zhang         /* perform dense axpy */
656e900a757SHong Zhang         valtmp = pA[j];
657a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
658e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
659a9d06231SHong Zhang         }
660a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
661a9d06231SHong Zhang       }
662a9d06231SHong Zhang 
663a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
664a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
665e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
666a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
667a9d06231SHong Zhang       for (j=0; j<pnz; j++){
668e900a757SHong Zhang         row = pdJ[j];
669e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
670e900a757SHong Zhang         cj  = bj + bi[row];
671e900a757SHong Zhang         ca  = ba + bi[row];
672a9d06231SHong Zhang         /* perform dense axpy */
673e900a757SHong Zhang         valtmp = pA[j];
674a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
675e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
676a9d06231SHong Zhang         }
677a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
678a9d06231SHong Zhang       }
679a9d06231SHong Zhang 
680a9d06231SHong Zhang       /* zero the current row of A*P */
681a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
682a9d06231SHong Zhang     }
683d9824c15SHong Zhang   } else {/* Perform sparse axpy */
684a9d06231SHong Zhang     /*----------------------------------------------------*/
685414894bdSHong Zhang     /* apa (length ap_rmax) stores sparse row A[i,:]*P */
686d6ab1dc0SHong Zhang     ierr = PetscMemzero(apa,ptap->rmax*sizeof(MatScalar));CHKERRQ(ierr);
68742c57489SHong Zhang 
688a9d06231SHong Zhang     pA=pa_loc;
68942c57489SHong Zhang     for (i=0; i<am; i++) {
690f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
69182412ba7SHong Zhang       apnz = api[i+1] - api[i];
69282412ba7SHong Zhang       apJ  = apj + api[i];
69342c57489SHong Zhang       /* diagonal portion of A */
69482412ba7SHong Zhang       anz = adi[i+1] - adi[i];
6951d633602SHong Zhang       adj = ad->j + adi[i];
6961d633602SHong Zhang       ada = ad->a + adi[i];
69782412ba7SHong Zhang       for (j=0; j<anz; j++) {
6981d633602SHong Zhang         row = adj[j];
69982412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
70082412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
70182412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
702e900a757SHong Zhang         valtmp = ada[j];
703f4a743e1SHong Zhang         nextp  = 0;
70482412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
70582412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
706e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
70742c57489SHong Zhang           }
70842c57489SHong Zhang         }
709dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
71042c57489SHong Zhang       }
71142c57489SHong Zhang       /* off-diagonal portion of A */
71282412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
7131d633602SHong Zhang       aoj = ao->j + aoi[i];
7141d633602SHong Zhang       aoa = ao->a + aoi[i];
71582412ba7SHong Zhang       for (j=0; j<anz; j++) {
7161d633602SHong Zhang         row = aoj[j];
71782412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
71882412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
71982412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
720e900a757SHong Zhang         valtmp = aoa[j];
721f4a743e1SHong Zhang         nextp  = 0;
72282412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
72382412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
724e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
72542c57489SHong Zhang           }
72642c57489SHong Zhang         }
727dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
72842c57489SHong Zhang       }
72942c57489SHong Zhang 
730a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
731a9d06231SHong Zhang       /*--------------------------------------------------------------*/
73282412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
733f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
73482412ba7SHong Zhang       for (j=0; j<pnz; j++) {
73542c57489SHong Zhang         nextap = 0;
736f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
73782412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
738e900a757SHong Zhang           row = *poJ;
739e900a757SHong Zhang           cj  = coj + coi[row];
740e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
74182412ba7SHong Zhang         } else {                            /* put the value into Cd */
742e900a757SHong Zhang           row  = *pdJ;
743e900a757SHong Zhang           cj   = bj + bi[row];
744e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
74542c57489SHong Zhang         }
746e900a757SHong Zhang         valtmp = pA[j];
74782412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
748e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
74942c57489SHong Zhang         }
750dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
751de0260b3SHong Zhang       }
7520496f32cSHong Zhang       pA += pnz;
75342c57489SHong Zhang       /* zero the current row info for A*P */
75482412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
75542c57489SHong Zhang     }
756d9824c15SHong Zhang   }
75742c57489SHong Zhang 
758a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
759a9d06231SHong Zhang   /*------------------------------------*/
76042c57489SHong Zhang   buf_ri = merge->buf_ri;
76142c57489SHong Zhang   buf_rj = merge->buf_rj;
76242c57489SHong Zhang   len_s  = merge->len_s;
763fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
76442c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
76542c57489SHong Zhang 
766a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
76742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
76842c57489SHong Zhang     if (!len_s[proc]) continue;
769de0260b3SHong Zhang     i = merge->owners_co[proc];
770de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
77142c57489SHong Zhang     k++;
77242c57489SHong Zhang   }
7730c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
7740c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
77542c57489SHong Zhang 
776a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
77742c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
77882412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
77942c57489SHong Zhang 
780a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
781a9d06231SHong Zhang   /*------------------------------------------------------*/
782687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
78342c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
78442c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
78542c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
78642c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
78782412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
78842c57489SHong Zhang   }
78942c57489SHong Zhang 
790de0260b3SHong Zhang   for (i=0; i<cm; i++) {
79182412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
79242c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
793de0260b3SHong Zhang     ba_i = ba + bi[i];
79482412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
79542c57489SHong Zhang     /* add received vals into ba */
79642c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
79742c57489SHong Zhang       /* i-th row */
79842c57489SHong Zhang       if (i == *nextrow[k]) {
79982412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
80082412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
80182412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
80282412ba7SHong Zhang         nextcj = 0;
80382412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
80482412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
80582412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
80642c57489SHong Zhang           }
80742c57489SHong Zhang         }
80882412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
809c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
81042c57489SHong Zhang       }
81142c57489SHong Zhang     }
81282412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
81342c57489SHong Zhang   }
81442c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81542c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81642c57489SHong Zhang 
817de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
818c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
81942c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
8200572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
821d9824c15SHong Zhang #if defined(PETSC_USE_INFO)
822d9824c15SHong Zhang   if (scalable){
823d9824c15SHong Zhang     ierr = PetscInfo(C,"Use scalable sparse axpy\n");CHKERRQ(ierr);
824d9824c15SHong Zhang   } else {
825d9824c15SHong Zhang     ierr = PetscInfo(C,"Use non-scalable dense axpy\n");CHKERRQ(ierr);
826d9824c15SHong Zhang   }
827d9824c15SHong Zhang #endif
82842c57489SHong Zhang   PetscFunctionReturn(0);
82942c57489SHong Zhang }
830