xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 5cab4f04ef148687507f6caacff5c5d88946249c)
1be1d678aSKris Buschelman 
22499ec78SHong Zhang /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
42499ec78SHong Zhang           C = A * B
52499ec78SHong Zhang */
6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
8c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>
11af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
121634b032SHong Zhang 
135e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
145e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
155e5acdf2Sstefano_zampini #endif
165e5acdf2Sstefano_zampini 
17150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
182499ec78SHong Zhang {
192499ec78SHong Zhang   PetscErrorCode ierr;
205e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
215e5acdf2Sstefano_zampini   const char     *algTypes[3] = {"scalable","nonscalable","hypre"};
225e5acdf2Sstefano_zampini   PetscInt       nalg = 3;
235e5acdf2Sstefano_zampini #else
240fc8cf34SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
255e5acdf2Sstefano_zampini   PetscInt       nalg = 2;
265e5acdf2Sstefano_zampini #endif
270d3441aeSHong Zhang   PetscInt       alg = 1; /* set default algorithm */
28e55be12dSHong Zhang   MPI_Comm       comm;
292499ec78SHong Zhang 
302499ec78SHong Zhang   PetscFunctionBegin;
312499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
32e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
336c4ed002SBarry Smith     if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
34e55be12dSHong Zhang 
350fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
365e5acdf2Sstefano_zampini     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,NULL);CHKERRQ(ierr);
370fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
380fc8cf34SHong Zhang 
393ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
400fc8cf34SHong Zhang     switch (alg) {
410fc8cf34SHong Zhang     case 1:
420fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
430fc8cf34SHong Zhang       break;
445e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
455e5acdf2Sstefano_zampini     case 2:
465e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
475e5acdf2Sstefano_zampini       break;
485e5acdf2Sstefano_zampini #endif
490fc8cf34SHong Zhang     default:
50b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
510fc8cf34SHong Zhang       break;
520fc8cf34SHong Zhang     }
533ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
542499ec78SHong Zhang   }
553ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
56598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
573ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
582499ec78SHong Zhang   PetscFunctionReturn(0);
592499ec78SHong Zhang }
602499ec78SHong Zhang 
61a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
62598bc09dSHong Zhang {
63598bc09dSHong Zhang   PetscErrorCode ierr;
64598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
65598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
66598bc09dSHong Zhang 
67598bc09dSHong Zhang   PetscFunctionBegin;
68b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
69598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
70a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
71a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
72ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
73a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
74a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
75598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
76598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
77598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
78598bc09dSHong Zhang   PetscFunctionReturn(0);
79598bc09dSHong Zhang }
80598bc09dSHong Zhang 
81a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
824ae0847bSHong Zhang {
834ae0847bSHong Zhang   PetscErrorCode ierr;
844ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
854ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
864ae0847bSHong Zhang 
874ae0847bSHong Zhang   PetscFunctionBegin;
884ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
8926fbe8dcSKarl Rupp 
904ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
914ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
924ae0847bSHong Zhang   PetscFunctionReturn(0);
934ae0847bSHong Zhang }
944ae0847bSHong Zhang 
95f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
96598bc09dSHong Zhang {
97598bc09dSHong Zhang   PetscErrorCode ierr;
984ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
99598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
100598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
101c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
102598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
1039ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
104c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
105598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
106c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
10729825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
108598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1099467f45fSHong Zhang   MPI_Comm       comm;
1109467f45fSHong Zhang   PetscMPIInt    size;
111598bc09dSHong Zhang 
112598bc09dSHong Zhang   PetscFunctionBegin;
1139467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1149467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1159467f45fSHong Zhang 
116a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
117598bc09dSHong Zhang   /*-----------------------------------------------------*/
118a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
119b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
120a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
121598bc09dSHong Zhang 
122598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
123598bc09dSHong Zhang   /*----------------------------------------------------------*/
124598bc09dSHong Zhang   /* get data from symbolic products */
125a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
126c12557a1SHong Zhang   p_oth = NULL;
1279467f45fSHong Zhang   if (size >1) {
1289467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1299467f45fSHong Zhang   }
130598bc09dSHong Zhang 
131598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
132598bc09dSHong Zhang   apa = ptap->apa;
133598bc09dSHong Zhang 
13429825010SHong Zhang   api = ptap->api;
13529825010SHong Zhang   apj = ptap->apj;
136598bc09dSHong Zhang   for (i=0; i<cm; i++) {
137c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
138e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
139598bc09dSHong Zhang 
140598bc09dSHong Zhang     /* set values in C */
141598bc09dSHong Zhang     apJ  = apj + api[i];
142598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
143598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
144598bc09dSHong Zhang 
145598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
146598bc09dSHong Zhang     ca = coa + co->i[i];
147598bc09dSHong Zhang     k  = 0;
148598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
149598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
150598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
151*5cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
152598bc09dSHong Zhang     }
153598bc09dSHong Zhang 
154598bc09dSHong Zhang     /* diagonal part of C */
155598bc09dSHong Zhang     ca = cda + cd->i[i];
156598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
157598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
158*5cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
159598bc09dSHong Zhang     }
160598bc09dSHong Zhang 
161598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
162598bc09dSHong Zhang     ca = coa + co->i[i];
163598bc09dSHong Zhang     for (; k0<conz; k0++) {
164598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
165*5cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
166598bc09dSHong Zhang     }
167598bc09dSHong Zhang   }
168598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170598bc09dSHong Zhang   PetscFunctionReturn(0);
171598bc09dSHong Zhang }
172598bc09dSHong Zhang 
1730fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
174598bc09dSHong Zhang {
175598bc09dSHong Zhang   PetscErrorCode     ierr;
176ce94432eSBarry Smith   MPI_Comm           comm;
1779467f45fSHong Zhang   PetscMPIInt        size;
178bfcd1a73SHong Zhang   Mat                Cmpi;
179598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1800298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1814ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
182bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1834ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1844ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
185d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
186*5cab4f04SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
187598bc09dSHong Zhang   PetscBT            lnkbt;
188598bc09dSHong Zhang   PetscScalar        *apa;
189bfcd1a73SHong Zhang   PetscReal          afill;
190598bc09dSHong Zhang 
191598bc09dSHong Zhang   PetscFunctionBegin;
192ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1939467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1949467f45fSHong Zhang 
195598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
196b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
197598bc09dSHong Zhang 
198598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
199b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
20019f0e57aSHong Zhang 
201598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
202a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
203598bc09dSHong Zhang 
204a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
205598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2069467f45fSHong Zhang   if (size > 1) {
2079467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
208598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
20920e3dc75SHong Zhang   } else {
21020e3dc75SHong Zhang     p_oth = NULL;
21120e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2129467f45fSHong Zhang   }
213598bc09dSHong Zhang 
214598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
215598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
216854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
217a1a4e84aSHong Zhang   ptap->api = api;
218598bc09dSHong Zhang   api[0]    = 0;
219598bc09dSHong Zhang 
220*5cab4f04SHong Zhang   /* create and initialize a linked list */
221*5cab4f04SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
222598bc09dSHong Zhang 
223bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
224f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
225598bc09dSHong Zhang   current_space = free_space;
226598bc09dSHong Zhang 
227598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
228598bc09dSHong Zhang   for (i=0; i<am; i++) {
229598bc09dSHong Zhang     /* diagonal portion of A */
230598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
231598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
232598bc09dSHong Zhang       row  = *adj++;
233598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
234598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
235598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2361fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
237598bc09dSHong Zhang     }
238598bc09dSHong Zhang     /* off-diagonal portion of A */
239598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
240598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
241598bc09dSHong Zhang       row  = *aoj++;
242598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
243598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2441fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
245598bc09dSHong Zhang     }
246598bc09dSHong Zhang 
247d14fa243SHong Zhang     apnz     = lnk[0];
248598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
249598bc09dSHong Zhang 
250598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
251598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
252f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
253598bc09dSHong Zhang       nspacedouble++;
254598bc09dSHong Zhang     }
255598bc09dSHong Zhang 
256598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2571fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
258598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2592205254eSKarl Rupp 
260598bc09dSHong Zhang     current_space->array           += apnz;
261598bc09dSHong Zhang     current_space->local_used      += apnz;
262598bc09dSHong Zhang     current_space->local_remaining -= apnz;
263598bc09dSHong Zhang   }
264598bc09dSHong Zhang 
265598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
266598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
267854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
268a1a4e84aSHong Zhang   apj  = ptap->apj;
269a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
270598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
271598bc09dSHong Zhang 
272f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2731795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2742205254eSKarl Rupp 
275f84c536bSHong Zhang   ptap->apa = apa;
276f84c536bSHong Zhang 
277bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
278598bc09dSHong Zhang   /*----------------------------------------------------*/
279bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
280bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
28133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
282a2f3521dSMark F. Adams 
283bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
284bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
285598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
286598bc09dSHong Zhang   for (i=0; i<am; i++) {
287598bc09dSHong Zhang     row  = i + rstart;
288598bc09dSHong Zhang     apnz = api[i+1] - api[i];
289bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
290598bc09dSHong Zhang     apj += apnz;
291598bc09dSHong Zhang   }
292bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
294598bc09dSHong Zhang 
295bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
296bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
297f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
298bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
299bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
300598bc09dSHong Zhang 
301bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
302bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
303598bc09dSHong Zhang   c->ptap = ptap;
304598bc09dSHong Zhang 
305bfcd1a73SHong Zhang   *C = Cmpi;
306bfcd1a73SHong Zhang 
307bfcd1a73SHong Zhang   /* set MatInfo */
308118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
309bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
310bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
311bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
312bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
313bfcd1a73SHong Zhang 
314bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
315bfcd1a73SHong Zhang   if (api[am]) {
31657622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
31757622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
318bfcd1a73SHong Zhang   } else {
319bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
320bfcd1a73SHong Zhang   }
321bfcd1a73SHong Zhang #endif
322598bc09dSHong Zhang   PetscFunctionReturn(0);
323598bc09dSHong Zhang }
324598bc09dSHong Zhang 
325150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3269123193aSHong Zhang {
3279123193aSHong Zhang   PetscErrorCode ierr;
3289123193aSHong Zhang 
3299123193aSHong Zhang   PetscFunctionBegin;
3309123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3313ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3329123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3333ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3349123193aSHong Zhang   }
3353ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3369123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3373ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3389123193aSHong Zhang   PetscFunctionReturn(0);
3399123193aSHong Zhang }
3409123193aSHong Zhang 
3414324174eSBarry Smith typedef struct {
3424324174eSBarry Smith   Mat         workB;
3434324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3444324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3454324174eSBarry Smith } MPIAIJ_MPIDense;
3464324174eSBarry Smith 
34796b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3484324174eSBarry Smith {
3494324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3504324174eSBarry Smith   PetscErrorCode  ierr;
3514324174eSBarry Smith 
3524324174eSBarry Smith   PetscFunctionBegin;
3536bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3544324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3554324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3564324174eSBarry Smith   PetscFunctionReturn(0);
3574324174eSBarry Smith }
3584324174eSBarry Smith 
359fd4e9aacSBarry Smith /*
360fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
361fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
362fd4e9aacSBarry Smith 
363fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
364fd4e9aacSBarry Smith */
365fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
366fd4e9aacSBarry Smith {
367fd4e9aacSBarry Smith   PetscErrorCode         ierr;
368fd4e9aacSBarry Smith   PetscBool              flg;
369fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
370fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
371fd4e9aacSBarry Smith   PetscContainer         container;
372fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
373fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
374fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
375fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
376fd4e9aacSBarry Smith 
377fd4e9aacSBarry Smith   PetscFunctionBegin;
378fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
379fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
380fd4e9aacSBarry Smith 
381fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
382fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
383fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
384fd4e9aacSBarry Smith 
385fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
386fd4e9aacSBarry Smith 
387fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
388fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
389fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
390fd4e9aacSBarry Smith   /* Create work arrays needed */
391fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
392fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
393fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
394fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
395fd4e9aacSBarry Smith 
396fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
397fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
398fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
399fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
400fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
401fd4e9aacSBarry Smith 
402fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
403fd4e9aacSBarry Smith   PetscFunctionReturn(0);
404fd4e9aacSBarry Smith }
405fd4e9aacSBarry Smith 
4069123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4079123193aSHong Zhang {
4089123193aSHong Zhang   PetscErrorCode         ierr;
409f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
410d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
411bf0cc555SLisandro Dalcin   PetscContainer         container;
4124324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4134324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4144324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4154324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
416d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4179123193aSHong Zhang 
4189123193aSHong Zhang   PetscFunctionBegin;
419ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
420d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
42133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
422cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4230298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
424cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
425cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4262205254eSKarl Rupp 
427d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
428f32d5d43SBarry Smith 
429b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
430f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4310298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4324324174eSBarry Smith   /* Create work arrays needed */
433dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
434dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
435dcca6d9dSJed Brown                       from->n,&contents->rwaits,
436dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4374324174eSBarry Smith 
438ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
439bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
44096b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
441bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
442bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4439123193aSHong Zhang   PetscFunctionReturn(0);
4449123193aSHong Zhang }
4459123193aSHong Zhang 
446f32d5d43SBarry Smith /*
4472636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4482636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
449f32d5d43SBarry Smith */
4504324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
451f32d5d43SBarry Smith {
452f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
453f32d5d43SBarry Smith   PetscErrorCode         ierr;
454f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
455f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
456f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
457f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
458f32d5d43SBarry Smith   PetscInt               i,j,k;
459aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
460aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
461f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
462ce94432eSBarry Smith   MPI_Comm               comm;
463d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
464f32d5d43SBarry Smith   MPI_Status             status;
4654324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
466bf0cc555SLisandro Dalcin   PetscContainer         container;
4674324174eSBarry Smith   Mat                    workB;
468f32d5d43SBarry Smith 
469f32d5d43SBarry Smith   PetscFunctionBegin;
470ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
471bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
47229825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
473bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4744324174eSBarry Smith 
4754324174eSBarry Smith   workB = *outworkB = contents->workB;
476cf1a0672SHong Zhang   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
477f32d5d43SBarry Smith   sindices = to->indices;
478f32d5d43SBarry Smith   sstarts  = to->starts;
479f32d5d43SBarry Smith   sprocs   = to->procs;
4804324174eSBarry Smith   swaits   = contents->swaits;
4814324174eSBarry Smith   svalues  = contents->svalues;
482f32d5d43SBarry Smith 
483f32d5d43SBarry Smith   rindices = from->indices;
484f32d5d43SBarry Smith   rstarts  = from->starts;
485f32d5d43SBarry Smith   rprocs   = from->procs;
4864324174eSBarry Smith   rwaits   = contents->rwaits;
4874324174eSBarry Smith   rvalues  = contents->rvalues;
488f32d5d43SBarry Smith 
4898c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
4908c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
491f32d5d43SBarry Smith 
492f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
493f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
494f32d5d43SBarry Smith   }
495f32d5d43SBarry Smith 
496f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
497f32d5d43SBarry Smith     /* pack a message at a time */
498f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
499f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
500f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5012636ff26SBarry Smith       }
5022636ff26SBarry Smith     }
503f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
504f32d5d43SBarry Smith   }
5052636ff26SBarry Smith 
506f32d5d43SBarry Smith   nrecvs = from->n;
507f32d5d43SBarry Smith   while (nrecvs) {
508f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
509f32d5d43SBarry Smith     nrecvs--;
510f32d5d43SBarry Smith     /* unpack a message at a time */
511f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
512f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
513f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5142636ff26SBarry Smith       }
5152636ff26SBarry Smith     }
516f32d5d43SBarry Smith   }
517cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
518f32d5d43SBarry Smith 
5198c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5208c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
521f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
522f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
523f32d5d43SBarry Smith   PetscFunctionReturn(0);
524f32d5d43SBarry Smith }
525f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
526f32d5d43SBarry Smith 
5279123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5289123193aSHong Zhang {
5299123193aSHong Zhang   PetscErrorCode ierr;
530f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
531f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
532f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
533f32d5d43SBarry Smith   Mat            workB;
5349123193aSHong Zhang 
5359123193aSHong Zhang   PetscFunctionBegin;
536f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
537f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
538f32d5d43SBarry Smith 
539f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5404324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
541f32d5d43SBarry Smith 
542f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
543f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5449123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5459123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5469123193aSHong Zhang   PetscFunctionReturn(0);
5479123193aSHong Zhang }
548cf1a0672SHong Zhang 
549f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5501634b032SHong Zhang {
5511634b032SHong Zhang   PetscErrorCode ierr;
552cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
553cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
554cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
555cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
556cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
557cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
558cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
55929825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
56029825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
561cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
56229825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
563cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
56429825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
56529825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
566a7c7454dSHong Zhang   MPI_Comm       comm;
567a7c7454dSHong Zhang   PetscMPIInt    size;
5681634b032SHong Zhang 
5691634b032SHong Zhang   PetscFunctionBegin;
570a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
571a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
572a7c7454dSHong Zhang 
573cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
574cf1a0672SHong Zhang   /*-----------------------------------------------------*/
575cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
576b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
577cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
578cf1a0672SHong Zhang 
579cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
580cf1a0672SHong Zhang   /*----------------------------------------------------------*/
581cf1a0672SHong Zhang   /* get data from symbolic products */
582cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
583cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
584a7c7454dSHong Zhang   if (size >1) {
585a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
586cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
58720e3dc75SHong Zhang   } else {
58820e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
589a7c7454dSHong Zhang   }
590cf1a0672SHong Zhang 
59129825010SHong Zhang   api = ptap->api;
59229825010SHong Zhang   apj = ptap->apj;
593cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
59429825010SHong Zhang     apJ = apj + api[i];
59529825010SHong Zhang 
596cf1a0672SHong Zhang     /* diagonal portion of A */
597cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
598cf1a0672SHong Zhang     adj = ad->j + adi[i];
599cf1a0672SHong Zhang     ada = ad->a + adi[i];
600cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
601cf1a0672SHong Zhang       row = adj[j];
602cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
603cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
604cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
60529825010SHong Zhang       /* perform sparse axpy */
606cf1a0672SHong Zhang       valtmp = ada[j];
60729825010SHong Zhang       nextp  = 0;
60829825010SHong Zhang       for (k=0; nextp<pnz; k++) {
60929825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
61029825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
61129825010SHong Zhang         }
6121634b032SHong Zhang       }
613cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
614cf1a0672SHong Zhang     }
6151634b032SHong Zhang 
616cf1a0672SHong Zhang     /* off-diagonal portion of A */
617cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
618cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
619cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
620cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
621cf1a0672SHong Zhang       row = aoj[j];
622cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
623cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
624cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
62529825010SHong Zhang       /* perform sparse axpy */
626cf1a0672SHong Zhang       valtmp = aoa[j];
62729825010SHong Zhang       nextp  = 0;
62829825010SHong Zhang       for (k=0; nextp<pnz; k++) {
62929825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
63029825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63129825010SHong Zhang         }
632cf1a0672SHong Zhang       }
633cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
634cf1a0672SHong Zhang     }
6351634b032SHong Zhang 
636cf1a0672SHong Zhang     /* set values in C */
637cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
638cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6391634b032SHong Zhang 
640cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
641cf1a0672SHong Zhang     ca = coa + co->i[i];
642cf1a0672SHong Zhang     k  = 0;
643cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
644cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
64529825010SHong Zhang       ca[k0]        = apa_sparse[k];
64629825010SHong Zhang       apa_sparse[k] = 0.0;
647cf1a0672SHong Zhang       k++;
648cf1a0672SHong Zhang     }
6491634b032SHong Zhang 
650cf1a0672SHong Zhang     /* diagonal part of C */
651cf1a0672SHong Zhang     ca = cda + cd->i[i];
652cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
65329825010SHong Zhang       ca[k1]        = apa_sparse[k];
65429825010SHong Zhang       apa_sparse[k] = 0.0;
655cf1a0672SHong Zhang       k++;
656cf1a0672SHong Zhang     }
657cf1a0672SHong Zhang 
658cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
659cf1a0672SHong Zhang     ca = coa + co->i[i];
660cf1a0672SHong Zhang     for (; k0<conz; k0++) {
66129825010SHong Zhang       ca[k0]        = apa_sparse[k];
66229825010SHong Zhang       apa_sparse[k] = 0.0;
663cf1a0672SHong Zhang       k++;
664cf1a0672SHong Zhang     }
665cf1a0672SHong Zhang   }
666cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
667cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
668cf1a0672SHong Zhang   PetscFunctionReturn(0);
669cf1a0672SHong Zhang }
670cf1a0672SHong Zhang 
6710fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
672b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
673cf1a0672SHong Zhang {
674cf1a0672SHong Zhang   PetscErrorCode     ierr;
675ce94432eSBarry Smith   MPI_Comm           comm;
676a7c7454dSHong Zhang   PetscMPIInt        size;
677cf1a0672SHong Zhang   Mat                Cmpi;
678cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6790298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
680cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
681cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
682cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
683cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
6840ca7d551SHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max;
6850ca7d551SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
686cf1a0672SHong Zhang   PetscReal          afill;
687cf1a0672SHong Zhang   PetscScalar        *apa;
688eca6b86aSHong Zhang   PetscTable         ta;
689cf1a0672SHong Zhang 
690cf1a0672SHong Zhang   PetscFunctionBegin;
691ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
692a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
693a7c7454dSHong Zhang 
694cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
695b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
696cf1a0672SHong Zhang 
697cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
698b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
69919f0e57aSHong Zhang 
700cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
701cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
702cf1a0672SHong Zhang 
703cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
704cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
705a7c7454dSHong Zhang   if (size > 1) {
706a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
70720e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
708968382fdSHong Zhang   } else {
70920e3dc75SHong Zhang     p_oth  = NULL;
71020e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
711a7c7454dSHong Zhang   }
712cf1a0672SHong Zhang 
713cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
714cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
715854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
716cf1a0672SHong Zhang   ptap->api = api;
717cf1a0672SHong Zhang   api[0]    = 0;
718cf1a0672SHong Zhang 
719cf1a0672SHong Zhang   /* create and initialize a linked list */
72045d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr);
7210ca7d551SHong Zhang 
7220ca7d551SHong Zhang   /* Calculate apnz_max */
7230ca7d551SHong Zhang   apnz_max = 0;
7240ca7d551SHong Zhang   for (i=0; i<am; i++) {
7250ca7d551SHong Zhang     ierr = PetscTableRemoveAll(ta);CHKERRQ(ierr);
7260ca7d551SHong Zhang     /* diagonal portion of A */
7270ca7d551SHong Zhang     nzi  = adi[i+1] - adi[i];
7280ca7d551SHong Zhang     Jptr = adj+adi[i];  /* cols of A_diag */
7290ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_loc,nzi,Jptr,ta);
7300ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7310ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7320ca7d551SHong Zhang 
7330ca7d551SHong Zhang     /*  off-diagonal portion of A */
7340ca7d551SHong Zhang     nzi = aoi[i+1] - aoi[i];
7350ca7d551SHong Zhang     Jptr = aoj+aoi[i];  /* cols of A_off */
7360ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_oth,nzi,Jptr,ta);
7370ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7380ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7390ca7d551SHong Zhang   }
740eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
741eca6b86aSHong Zhang 
7420ca7d551SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(apnz_max,&lnk);CHKERRQ(ierr);
743cf1a0672SHong Zhang 
744cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
745f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
746cf1a0672SHong Zhang   current_space = free_space;
747cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
748cf1a0672SHong Zhang   for (i=0; i<am; i++) {
749cf1a0672SHong Zhang     /* diagonal portion of A */
750cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
751cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
752cf1a0672SHong Zhang       row  = *adj++;
753cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
754cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
755cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
756f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
757cf1a0672SHong Zhang     }
758cf1a0672SHong Zhang     /* off-diagonal portion of A */
759cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
760cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
761cf1a0672SHong Zhang       row  = *aoj++;
762cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
763cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
764f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
765cf1a0672SHong Zhang     }
766cf1a0672SHong Zhang 
767f84c536bSHong Zhang     apnz     = *lnk;
768cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
769cf1a0672SHong Zhang 
770cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
771cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
772f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
773cf1a0672SHong Zhang       nspacedouble++;
774cf1a0672SHong Zhang     }
775cf1a0672SHong Zhang 
776cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
777f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
778cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7792205254eSKarl Rupp 
780cf1a0672SHong Zhang     current_space->array           += apnz;
781cf1a0672SHong Zhang     current_space->local_used      += apnz;
782cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
783cf1a0672SHong Zhang   }
784cf1a0672SHong Zhang 
785cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
786cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
787854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
788cf1a0672SHong Zhang   apj  = ptap->apj;
789cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
790f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
791cf1a0672SHong Zhang 
792cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
793cf1a0672SHong Zhang   /*----------------------------------------------------*/
794cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
795cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
79633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
797cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
798cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
799cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
800cf1a0672SHong Zhang 
80129825010SHong Zhang   /* malloc apa for assembly Cmpi */
8021795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8032205254eSKarl Rupp 
80429825010SHong Zhang   ptap->apa = apa;
805cf1a0672SHong Zhang   for (i=0; i<am; i++) {
806cf1a0672SHong Zhang     row  = i + rstart;
807cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
808cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
809cf1a0672SHong Zhang     apj += apnz;
810cf1a0672SHong Zhang   }
811cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
812cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
813cf1a0672SHong Zhang 
814cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
815cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
816f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
817cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
818cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
819cf1a0672SHong Zhang 
820cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
821cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
822cf1a0672SHong Zhang   c->ptap = ptap;
823cf1a0672SHong Zhang 
824cf1a0672SHong Zhang   *C = Cmpi;
825cf1a0672SHong Zhang 
826cf1a0672SHong Zhang   /* set MatInfo */
827118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
828cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
829cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
830cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
831cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
832cf1a0672SHong Zhang 
833cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
834cf1a0672SHong Zhang   if (api[am]) {
83557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
83657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
837cf1a0672SHong Zhang   } else {
838cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
839cf1a0672SHong Zhang   }
840cf1a0672SHong Zhang #endif
8411634b032SHong Zhang   PetscFunctionReturn(0);
8421634b032SHong Zhang }
843f96d37fcSHong Zhang 
844f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
845c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
846f96d37fcSHong Zhang {
847f96d37fcSHong Zhang   PetscErrorCode ierr;
848c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
849c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
850f96d37fcSHong Zhang 
851f96d37fcSHong Zhang   PetscFunctionBegin;
852c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
853c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
854c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
855c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
856c216dbf3SHong Zhang 
857c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
858c216dbf3SHong Zhang     switch (alg) {
859c216dbf3SHong Zhang     case 1:
8602bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
861c216dbf3SHong Zhang       break;
862c216dbf3SHong Zhang     case 2:
863275476c6SMatthew G. Knepley     {
864ab1b013aSHong Zhang       Mat         Pt;
865ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
866ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
867ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
868ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
869ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
870ab1b013aSHong Zhang       ptap     = c->ptap;
871ab1b013aSHong Zhang       ptap->Pt = Pt;
872c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
873c216dbf3SHong Zhang       PetscFunctionReturn(0);
874275476c6SMatthew G. Knepley     }
875c216dbf3SHong Zhang       break;
876c216dbf3SHong Zhang     default:
8776da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
878c216dbf3SHong Zhang       break;
879c216dbf3SHong Zhang     }
880c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
881c216dbf3SHong Zhang   }
882c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
883c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
884c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
885ab1b013aSHong Zhang   PetscFunctionReturn(0);
886ab1b013aSHong Zhang }
887ab1b013aSHong Zhang 
888c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
889c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
890c216dbf3SHong Zhang {
891c216dbf3SHong Zhang   PetscErrorCode ierr;
8922bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
8932bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
8942bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
895c216dbf3SHong Zhang 
896c216dbf3SHong Zhang   PetscFunctionBegin;
897c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
898c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
899f96d37fcSHong Zhang   PetscFunctionReturn(0);
900f96d37fcSHong Zhang }
901f96d37fcSHong Zhang 
9026da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
9036da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
904f96d37fcSHong Zhang {
905c5af039cSHong Zhang   PetscErrorCode      ierr;
906c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
907497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
908c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
909c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
910d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
911497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
912e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
913c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
914ce94432eSBarry Smith   MPI_Comm            comm;
915c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
916c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
917c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
918c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
919c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
920c5af039cSHong Zhang   MPI_Status          *status;
921c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
922d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
923a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
924c5af039cSHong Zhang   Mat                 A_loc;
925c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
926f96d37fcSHong Zhang 
927f96d37fcSHong Zhang   PetscFunctionBegin;
928ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
929c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
930c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
931c5af039cSHong Zhang 
932c5af039cSHong Zhang   ptap  = c->ptap;
933c5af039cSHong Zhang   merge = ptap->merge;
934c5af039cSHong Zhang 
935c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
936c5af039cSHong Zhang   /*--------------------------------------------------------------*/
937c5af039cSHong Zhang   /* get data from symbolic products */
938c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
939854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
940c5af039cSHong Zhang 
941c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
942c5af039cSHong Zhang   owners = merge->rowmap->range;
943854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
944c5af039cSHong Zhang 
945c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
946c5af039cSHong Zhang   A_loc = ptap->A_loc;
947c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
948c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
949d6ab1dc0SHong Zhang   ai    = a_loc->i;
950d6ab1dc0SHong Zhang   aj    = a_loc->j;
951c5af039cSHong Zhang 
952854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
953c5af039cSHong Zhang 
954c5af039cSHong Zhang   for (i=0; i<am; i++) {
955e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
956d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
957d6ab1dc0SHong Zhang     adj = aj + ai[i];
958d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
959c5af039cSHong Zhang     for (j=0; j<anz; j++) {
960e745005bSHong Zhang       aval[adj[j]] = ada[j];
961c5af039cSHong Zhang     }
962c5af039cSHong Zhang 
963c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
964c5af039cSHong Zhang     /*--------------------------------------------------------------*/
965c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
966c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
967c5af039cSHong Zhang     poJ = po->j + po->i[i];
968c5af039cSHong Zhang     pA  = po->a + po->i[i];
969c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
970c5af039cSHong Zhang       row = poJ[j];
971c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
972c5af039cSHong Zhang       cj  = coj + coi[row];
973c5af039cSHong Zhang       ca  = coa + coi[row];
974c5af039cSHong Zhang       /* perform dense axpy */
975c5af039cSHong Zhang       valtmp = pA[j];
976c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
977e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
978c5af039cSHong Zhang       }
979c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
980c5af039cSHong Zhang     }
981c5af039cSHong Zhang 
982c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
983c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
984c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
985c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
986c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
987c5af039cSHong Zhang       row = pdJ[j];
988c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
989c5af039cSHong Zhang       cj  = bj + bi[row];
990c5af039cSHong Zhang       ca  = ba + bi[row];
991c5af039cSHong Zhang       /* perform dense axpy */
992c5af039cSHong Zhang       valtmp = pA[j];
993c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
994e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
995c5af039cSHong Zhang       }
996c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
997c5af039cSHong Zhang     }
998c5af039cSHong Zhang 
999d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1000d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1001e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1002c5af039cSHong Zhang   }
1003c5af039cSHong Zhang 
1004c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1005c5af039cSHong Zhang   /*------------------------------------*/
1006c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1007c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1008c5af039cSHong Zhang   len_s  = merge->len_s;
1009c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1010c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1011c5af039cSHong Zhang 
1012dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1013c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1014c5af039cSHong Zhang     if (!len_s[proc]) continue;
1015c5af039cSHong Zhang     i    = merge->owners_co[proc];
1016c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1017c5af039cSHong Zhang     k++;
1018c5af039cSHong Zhang   }
1019c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1020c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1021c5af039cSHong Zhang 
1022c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1023c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1024c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1025c5af039cSHong Zhang 
1026c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1027c5af039cSHong Zhang   /*----------------------------------------------------*/
1028dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1029c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1030c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1031c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1032c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1033c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1034c5af039cSHong Zhang   }
1035c5af039cSHong Zhang 
1036c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1037c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1038c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1039c5af039cSHong Zhang     ba_i = ba + bi[i];
1040c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1041c5af039cSHong Zhang     /* add received vals into ba */
1042c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1043c5af039cSHong Zhang       /* i-th row */
1044c5af039cSHong Zhang       if (i == *nextrow[k]) {
1045c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1046c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1047c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1048c5af039cSHong Zhang         nextcj = 0;
1049c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1050c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1051c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1052c5af039cSHong Zhang           }
1053c5af039cSHong Zhang         }
1054c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1055c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1056c5af039cSHong Zhang       }
1057c5af039cSHong Zhang     }
1058c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1059c5af039cSHong Zhang   }
1060c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1061c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1062c5af039cSHong Zhang 
1063c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1064c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1065c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1066c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1067e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1068f96d37fcSHong Zhang   PetscFunctionReturn(0);
1069f96d37fcSHong Zhang }
1070f96d37fcSHong Zhang 
1071c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1072c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
10732bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1074f96d37fcSHong Zhang {
1075f96d37fcSHong Zhang   PetscErrorCode      ierr;
10764e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1077f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10780298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1079f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1080f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1081f96d37fcSHong Zhang   PetscInt            nnz;
10824e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1083497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1084f96d37fcSHong Zhang   PetscBT             lnkbt;
1085ce94432eSBarry Smith   MPI_Comm            comm;
1086f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1087f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1088f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1089f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1090f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1091f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1092f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1093f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1094d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1095f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1096438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1097f96d37fcSHong Zhang   PetscScalar         *vals;
10984e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1099f96d37fcSHong Zhang 
1100f96d37fcSHong Zhang   PetscFunctionBegin;
1101ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1102c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
11036c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1104c5af039cSHong Zhang 
1105f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1106f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1107f96d37fcSHong Zhang 
1108f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1109b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1110f96d37fcSHong Zhang 
1111f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1112f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11132205254eSKarl Rupp 
1114c5af039cSHong Zhang   ptap->A_loc = A_loc;
11152205254eSKarl Rupp 
11161c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1117d6ab1dc0SHong Zhang   ai    = a_loc->i;
1118d6ab1dc0SHong Zhang   aj    = a_loc->j;
1119f96d37fcSHong Zhang 
1120f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1121f96d37fcSHong Zhang   /*----------------------------------------------------*/
11224e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11234e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11244e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11254e938277SHong Zhang 
11264e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11274e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11284e938277SHong Zhang   poti = pot->i; potj = pot->j;
1129f96d37fcSHong Zhang 
1130f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11312205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1132854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1133f96d37fcSHong Zhang   coi[0] = 0;
1134f96d37fcSHong Zhang 
1135f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1136f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1137a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1138f96d37fcSHong Zhang   current_space = free_space;
113919f0e57aSHong Zhang 
1140f96d37fcSHong Zhang   /* create and initialize a linked list */
1141438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1142f96d37fcSHong Zhang 
1143f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1144f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1145f96d37fcSHong Zhang     ptJ = potj + poti[i];
1146f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1147f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1148d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1149d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1150f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1151d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1152f96d37fcSHong Zhang     }
11534e938277SHong Zhang     nnz = lnk[0];
1154f96d37fcSHong Zhang 
1155f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1156f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1157f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1158f96d37fcSHong Zhang       nspacedouble++;
1159f96d37fcSHong Zhang     }
1160f96d37fcSHong Zhang 
1161f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11624e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11632205254eSKarl Rupp 
1164f96d37fcSHong Zhang     current_space->array           += nnz;
1165f96d37fcSHong Zhang     current_space->local_used      += nnz;
1166f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11672205254eSKarl Rupp 
1168f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1169f96d37fcSHong Zhang   }
1170f96d37fcSHong Zhang 
1171854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1172f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11732205254eSKarl Rupp 
1174118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1175f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1176f96d37fcSHong Zhang 
1177f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1178f96d37fcSHong Zhang   /*----------------------------------------------*/
1179f96d37fcSHong Zhang   /* determine row ownership */
1180b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1181f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11822205254eSKarl Rupp 
1183f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1184f96d37fcSHong Zhang   merge->rowmap->bs = 1;
11852205254eSKarl Rupp 
1186f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1187f96d37fcSHong Zhang   owners = merge->rowmap->range;
1188f96d37fcSHong Zhang 
1189f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
11901795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1191785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
11922205254eSKarl Rupp 
1193f96d37fcSHong Zhang   len_s        = merge->len_s;
1194f96d37fcSHong Zhang   merge->nsend = 0;
1195f96d37fcSHong Zhang 
1196854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1197f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1198f96d37fcSHong Zhang 
1199f96d37fcSHong Zhang   proc = 0;
1200f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1201f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1202f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1203f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1204f96d37fcSHong Zhang   }
1205f96d37fcSHong Zhang 
1206f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1207f96d37fcSHong Zhang   owners_co[0] = 0;
1208f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1209f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1210f96d37fcSHong Zhang     if (len_si[proc]) {
1211f96d37fcSHong Zhang       merge->nsend++;
1212f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1213f96d37fcSHong Zhang       len         += len_si[proc];
1214f96d37fcSHong Zhang     }
1215f96d37fcSHong Zhang   }
1216f96d37fcSHong Zhang 
1217f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12180298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1219f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1220f96d37fcSHong Zhang 
1221f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1222f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1223f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1224854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1225f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1226f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1227f96d37fcSHong Zhang     i    = owners_co[proc];
1228f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1229f96d37fcSHong Zhang     k++;
1230f96d37fcSHong Zhang   }
1231f96d37fcSHong Zhang 
1232f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1233785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1234f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1235c280ed6aSJed Brown     PetscMPIInt icompleted;
1236c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1237f96d37fcSHong Zhang   }
1238f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1239f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1240f96d37fcSHong Zhang 
1241f96d37fcSHong Zhang   /* send and recv coi */
1242f96d37fcSHong Zhang   /*-------------------*/
1243f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1244f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1245854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1246f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1247f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1248f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1249f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1250f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1251f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1252f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1253f96d37fcSHong Zhang     */
1254f96d37fcSHong Zhang     /*-------------------------------------------*/
1255f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1256f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1257f96d37fcSHong Zhang     buf_si[0]   = nrows;
1258f96d37fcSHong Zhang     buf_si_i[0] = 0;
1259f96d37fcSHong Zhang     nrows       = 0;
1260f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1261f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1262f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1263f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1264f96d37fcSHong Zhang       nrows++;
1265f96d37fcSHong Zhang     }
1266f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1267f96d37fcSHong Zhang     k++;
1268f96d37fcSHong Zhang     buf_si += len_si[proc];
1269f96d37fcSHong Zhang   }
1270f96d37fcSHong Zhang   i = merge->nrecv;
1271f96d37fcSHong Zhang   while (i--) {
1272c280ed6aSJed Brown     PetscMPIInt icompleted;
1273c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1274f96d37fcSHong Zhang   }
1275f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1276f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1277f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1278f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1279f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1280f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1281f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1282f96d37fcSHong Zhang 
1283f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1284f96d37fcSHong Zhang   /*------------------------------------------*/
1285f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1286854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1287f96d37fcSHong Zhang   bi[0] = 0;
1288f96d37fcSHong Zhang 
1289c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1290f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1291a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1292f96d37fcSHong Zhang   current_space = free_space;
1293f96d37fcSHong Zhang 
1294dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1295f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1296f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1297f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1298f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1299f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1300f96d37fcSHong Zhang   }
1301f96d37fcSHong Zhang 
13021c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1303f96d37fcSHong Zhang   rmax = 0;
1304f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1305f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1306f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1307f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1308f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1309f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1310d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1311d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1312f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1313d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1314f96d37fcSHong Zhang     }
13154e938277SHong Zhang 
1316f96d37fcSHong Zhang     /* add received col data into lnk */
1317f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1318f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1319f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1320f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13214e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1322f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1323f96d37fcSHong Zhang       }
1324f96d37fcSHong Zhang     }
13254e938277SHong Zhang     nnz = lnk[0];
1326f96d37fcSHong Zhang 
1327f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1328f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1329f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1330f96d37fcSHong Zhang       nspacedouble++;
1331f96d37fcSHong Zhang     }
1332f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13334e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1334f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13352205254eSKarl Rupp 
1336f96d37fcSHong Zhang     current_space->array           += nnz;
1337f96d37fcSHong Zhang     current_space->local_used      += nnz;
1338f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13392205254eSKarl Rupp 
1340f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1341f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1342f96d37fcSHong Zhang   }
1343f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1344f96d37fcSHong Zhang 
1345854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1346f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13472205254eSKarl Rupp 
1348118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1349f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1350d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13514e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13524e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1353f96d37fcSHong Zhang 
13541c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13551c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13561795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1357f96d37fcSHong Zhang 
1358f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1359f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
136033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1361f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1362f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1363f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1364f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1365f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1366f96d37fcSHong Zhang     row  = i + rstart;
1367f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1368f96d37fcSHong Zhang     Jptr = bj + bi[i];
1369f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1370f96d37fcSHong Zhang   }
1371f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1372f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1373f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1374f96d37fcSHong Zhang 
1375f96d37fcSHong Zhang   merge->bi        = bi;
1376f96d37fcSHong Zhang   merge->bj        = bj;
1377f96d37fcSHong Zhang   merge->coi       = coi;
1378f96d37fcSHong Zhang   merge->coj       = coj;
1379f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1380f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1381f96d37fcSHong Zhang   merge->owners_co = owners_co;
1382f96d37fcSHong Zhang 
1383f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1384f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1385f96d37fcSHong Zhang   c->ptap     = ptap;
13860298fd71SBarry Smith   ptap->api   = NULL;
13870298fd71SBarry Smith   ptap->apj   = NULL;
1388f96d37fcSHong Zhang   ptap->merge = merge;
1389375ed354SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1390375ed354SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1391375ed354SHong Zhang 
1392375ed354SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1393375ed354SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1394375ed354SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1395d6ab1dc0SHong Zhang 
1396d6ab1dc0SHong Zhang   *C = Cmpi;
1397d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1398d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
139957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
140057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1401d6ab1dc0SHong Zhang   } else {
1402d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1403d6ab1dc0SHong Zhang   }
1404d6ab1dc0SHong Zhang #endif
1405d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1406d6ab1dc0SHong Zhang }
1407d6ab1dc0SHong Zhang 
14086da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1409d6ab1dc0SHong Zhang {
1410d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1411d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1412d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1413d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1414d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1415e745005bSHong Zhang   PetscInt            *adj;
1416e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1417e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1418d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1419ce94432eSBarry Smith   MPI_Comm            comm;
1420d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1421d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1422d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1423d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1424d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1425d6ab1dc0SHong Zhang   MPI_Status          *status;
1426d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1427d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1428a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1429d6ab1dc0SHong Zhang   Mat                 A_loc;
1430d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1431d6ab1dc0SHong Zhang 
1432d6ab1dc0SHong Zhang   PetscFunctionBegin;
1433ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1434d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1435d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1436d6ab1dc0SHong Zhang 
1437d6ab1dc0SHong Zhang   ptap  = c->ptap;
1438d6ab1dc0SHong Zhang   merge = ptap->merge;
1439d6ab1dc0SHong Zhang 
1440e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1441e745005bSHong Zhang   /*------------------------------------------*/
1442d6ab1dc0SHong Zhang   /* get data from symbolic products */
1443d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1444854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1445d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1446d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14471795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1448d6ab1dc0SHong Zhang 
1449d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1450d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1451d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1452d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1453d6ab1dc0SHong Zhang   ai    = a_loc->i;
1454d6ab1dc0SHong Zhang   aj    = a_loc->j;
1455d6ab1dc0SHong Zhang 
1456d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1457d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1458d6ab1dc0SHong Zhang     adj = aj + ai[i];
1459d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1460d6ab1dc0SHong Zhang 
1461d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1462e745005bSHong Zhang     /*-------------------------------------------------------------*/
1463d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1464d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1465d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1466d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1467d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1468d6ab1dc0SHong Zhang       row = poJ[j];
1469d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1470d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1471e745005bSHong Zhang       /* perform sparse axpy */
1472e745005bSHong Zhang       nexta  = 0;
1473d6ab1dc0SHong Zhang       valtmp = pA[j];
1474e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1475e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1476e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1477e745005bSHong Zhang           nexta++;
1478d6ab1dc0SHong Zhang         }
1479e745005bSHong Zhang       }
1480e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1481d6ab1dc0SHong Zhang     }
1482d6ab1dc0SHong Zhang 
1483d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1484d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1485d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1486d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1487d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1488d6ab1dc0SHong Zhang       row = pdJ[j];
1489d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1490d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1491e745005bSHong Zhang       /* perform sparse axpy */
1492e745005bSHong Zhang       nexta  = 0;
1493d6ab1dc0SHong Zhang       valtmp = pA[j];
1494e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1495e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1496e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1497e745005bSHong Zhang           nexta++;
1498d6ab1dc0SHong Zhang         }
1499e745005bSHong Zhang       }
1500e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1501d6ab1dc0SHong Zhang     }
1502d6ab1dc0SHong Zhang   }
1503d6ab1dc0SHong Zhang 
1504d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1505d6ab1dc0SHong Zhang   /*------------------------------------*/
1506d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1507d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1508d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1509d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1510d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1511d6ab1dc0SHong Zhang 
1512dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1513d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1514d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1515d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1516d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1517d6ab1dc0SHong Zhang     k++;
1518d6ab1dc0SHong Zhang   }
1519d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1520d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1521d6ab1dc0SHong Zhang 
1522d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1523d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1524d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1525d6ab1dc0SHong Zhang 
1526d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1527d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1528dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1530e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1531d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1532d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1533d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1534d6ab1dc0SHong Zhang   }
1535d6ab1dc0SHong Zhang 
1536d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1537d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1538d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1539d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1540d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1541d6ab1dc0SHong Zhang     /* add received vals into ba */
1542d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1543d6ab1dc0SHong Zhang       /* i-th row */
1544d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1545d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1546d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1547d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1548d6ab1dc0SHong Zhang         nextcj = 0;
1549d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1550d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1551d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1552d6ab1dc0SHong Zhang           }
1553d6ab1dc0SHong Zhang         }
1554d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1555d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1556d6ab1dc0SHong Zhang       }
1557d6ab1dc0SHong Zhang     }
1558d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1559d6ab1dc0SHong Zhang   }
1560d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1561d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1562d6ab1dc0SHong Zhang 
1563d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1564d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1565d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1566d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1567d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1568d6ab1dc0SHong Zhang }
1569d6ab1dc0SHong Zhang 
1570c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
15712bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15722bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
15736da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1574d6ab1dc0SHong Zhang {
1575d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1576d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1577d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15780298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
157933ba5e2fSHong Zhang   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1580d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1581d6ab1dc0SHong Zhang   PetscInt            nnz;
1582d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1583d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1584ce94432eSBarry Smith   MPI_Comm            comm;
1585d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1586d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1587d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1588d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1589d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1590d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1591d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1592d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1593d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1594d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
15954b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1596d6ab1dc0SHong Zhang   PetscScalar         *vals;
1597d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
15980acc65b4SHong Zhang   PetscTable          ta;
1599d6ab1dc0SHong Zhang 
1600d6ab1dc0SHong Zhang   PetscFunctionBegin;
1601ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1602d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
16036c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1604d6ab1dc0SHong Zhang 
1605d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1606d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1607d6ab1dc0SHong Zhang 
1608d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1609b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1610d6ab1dc0SHong Zhang 
1611d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1612d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16132205254eSKarl Rupp 
1614d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1615d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1616d6ab1dc0SHong Zhang   ai          = a_loc->i;
1617d6ab1dc0SHong Zhang   aj          = a_loc->j;
1618d6ab1dc0SHong Zhang 
1619d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1620d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1621d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1622d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1623d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1624d6ab1dc0SHong Zhang 
1625d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1626d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1627d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1628d6ab1dc0SHong Zhang 
1629d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1630d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1631d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1632854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1633d6ab1dc0SHong Zhang   coi[0] = 0;
1634d6ab1dc0SHong Zhang 
1635d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1636f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1637a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1638d6ab1dc0SHong Zhang   current_space = free_space;
163919f0e57aSHong Zhang 
1640d6ab1dc0SHong Zhang   /* create and initialize a linked list */
164133ba5e2fSHong Zhang   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
16424b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
16430acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
164433ba5e2fSHong Zhang 
16450acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1646d6ab1dc0SHong Zhang 
1647d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1648d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1649d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1650d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1651d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1652d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1653d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1654d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1655d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1656d6ab1dc0SHong Zhang     }
1657d6ab1dc0SHong Zhang     nnz = lnk[0];
1658d6ab1dc0SHong Zhang 
1659d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1660d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1661f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1662d6ab1dc0SHong Zhang       nspacedouble++;
1663d6ab1dc0SHong Zhang     }
1664d6ab1dc0SHong Zhang 
1665d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1666d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16672205254eSKarl Rupp 
1668d6ab1dc0SHong Zhang     current_space->array           += nnz;
1669d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1670d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16712205254eSKarl Rupp 
1672d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1673d6ab1dc0SHong Zhang   }
1674d6ab1dc0SHong Zhang 
1675854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1676d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16770acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
16782205254eSKarl Rupp 
1679118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1680d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1681d6ab1dc0SHong Zhang 
1682d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1683d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1684d6ab1dc0SHong Zhang   /* determine row ownership */
1685b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1686d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
16872205254eSKarl Rupp 
1688d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1689d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
16902205254eSKarl Rupp 
1691d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1692d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1693d6ab1dc0SHong Zhang 
1694d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
16951795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1696785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
16972205254eSKarl Rupp 
1698d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1699d6ab1dc0SHong Zhang   merge->nsend = 0;
1700d6ab1dc0SHong Zhang 
1701854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1702d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1703d6ab1dc0SHong Zhang 
1704d6ab1dc0SHong Zhang   proc = 0;
1705d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1706d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1707d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1708d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1709d6ab1dc0SHong Zhang   }
1710d6ab1dc0SHong Zhang 
1711d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1712d6ab1dc0SHong Zhang   owners_co[0] = 0;
1713d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1714d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1715d6ab1dc0SHong Zhang     if (len_si[proc]) {
1716d6ab1dc0SHong Zhang       merge->nsend++;
1717d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1718d6ab1dc0SHong Zhang       len         += len_si[proc];
1719d6ab1dc0SHong Zhang     }
1720d6ab1dc0SHong Zhang   }
1721d6ab1dc0SHong Zhang 
1722d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17230298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1724d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1725d6ab1dc0SHong Zhang 
1726d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1727d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1728d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1729854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1730d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1731d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1732d6ab1dc0SHong Zhang     i    = owners_co[proc];
1733d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1734d6ab1dc0SHong Zhang     k++;
1735d6ab1dc0SHong Zhang   }
1736d6ab1dc0SHong Zhang 
1737d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1738785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1739d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1740c280ed6aSJed Brown     PetscMPIInt icompleted;
1741c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1742d6ab1dc0SHong Zhang   }
1743d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1745d6ab1dc0SHong Zhang 
17460acc65b4SHong Zhang   /* add received column indices into table to update Armax */
174733ba5e2fSHong Zhang   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
17480acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
17490acc65b4SHong Zhang     Jptr = buf_rj[k];
17500acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
17510acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
17520acc65b4SHong Zhang     }
17530acc65b4SHong Zhang   }
17540acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
175533ba5e2fSHong Zhang   /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */
17560acc65b4SHong Zhang 
1757d6ab1dc0SHong Zhang   /* send and recv coi */
1758d6ab1dc0SHong Zhang   /*-------------------*/
1759d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1760d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1761854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1762d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1763d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1764d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1765d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1766d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1767d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1768d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1769d6ab1dc0SHong Zhang     */
1770d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1771d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1772d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1773d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1774d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1775d6ab1dc0SHong Zhang     nrows       = 0;
1776d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1777d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1778d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1779d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1780d6ab1dc0SHong Zhang       nrows++;
1781d6ab1dc0SHong Zhang     }
1782d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1783d6ab1dc0SHong Zhang     k++;
1784d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1785d6ab1dc0SHong Zhang   }
1786d6ab1dc0SHong Zhang   i = merge->nrecv;
1787d6ab1dc0SHong Zhang   while (i--) {
1788c280ed6aSJed Brown     PetscMPIInt icompleted;
1789c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1790d6ab1dc0SHong Zhang   }
1791d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1792d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1793d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1794d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1795d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1796d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1797d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1798d6ab1dc0SHong Zhang 
1799d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1800d6ab1dc0SHong Zhang   /*------------------------------------------*/
1801d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1802854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1803d6ab1dc0SHong Zhang   bi[0] = 0;
1804d6ab1dc0SHong Zhang 
1805d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1806f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1807a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1808d6ab1dc0SHong Zhang   current_space = free_space;
1809d6ab1dc0SHong Zhang 
1810dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1811d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1812d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1813d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1814d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18152205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1816d6ab1dc0SHong Zhang   }
1817d6ab1dc0SHong Zhang 
18180acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1820d6ab1dc0SHong Zhang   rmax = 0;
1821d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1822d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1823d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1824d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1825d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1826d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1827d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1828d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1829d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1830d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1831d6ab1dc0SHong Zhang     }
1832d6ab1dc0SHong Zhang 
1833d6ab1dc0SHong Zhang     /* add received col data into lnk */
1834d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1835d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1836d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1837d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1838d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1839d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1840d6ab1dc0SHong Zhang       }
1841d6ab1dc0SHong Zhang     }
1842d6ab1dc0SHong Zhang     nnz = lnk[0];
1843d6ab1dc0SHong Zhang 
1844d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1845d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1846f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1847d6ab1dc0SHong Zhang       nspacedouble++;
1848d6ab1dc0SHong Zhang     }
1849d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1850d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1851d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18522205254eSKarl Rupp 
1853d6ab1dc0SHong Zhang     current_space->array           += nnz;
1854d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1855d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18562205254eSKarl Rupp 
1857d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1858d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1859d6ab1dc0SHong Zhang   }
1860f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1861d6ab1dc0SHong Zhang 
1862854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1863d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1864118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1865d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1866d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
18670acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
18680acc65b4SHong Zhang 
1869d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1870d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1871d6ab1dc0SHong Zhang 
1872d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1873d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18741795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1875d6ab1dc0SHong Zhang 
1876d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1877d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
187833d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1879d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1881d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1882d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1883d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1884d6ab1dc0SHong Zhang     row  = i + rstart;
1885d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1886d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1887d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1888d6ab1dc0SHong Zhang   }
1889d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1890d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1891d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1892d6ab1dc0SHong Zhang 
1893d6ab1dc0SHong Zhang   merge->bi        = bi;
1894d6ab1dc0SHong Zhang   merge->bj        = bj;
1895d6ab1dc0SHong Zhang   merge->coi       = coi;
1896d6ab1dc0SHong Zhang   merge->coj       = coj;
1897d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1898d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1899d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1900d6ab1dc0SHong Zhang 
1901d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1902d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19032205254eSKarl Rupp 
1904d6ab1dc0SHong Zhang   c->ptap     = ptap;
19050298fd71SBarry Smith   ptap->api   = NULL;
19060298fd71SBarry Smith   ptap->apj   = NULL;
1907d6ab1dc0SHong Zhang   ptap->merge = merge;
19080298fd71SBarry Smith   ptap->apa   = NULL;
1909a560ef98SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1910a560ef98SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1911a560ef98SHong Zhang 
1912a560ef98SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1913a560ef98SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1914a560ef98SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1915f96d37fcSHong Zhang 
1916f96d37fcSHong Zhang   *C = Cmpi;
1917f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1918f96d37fcSHong Zhang   if (bi[pn] != 0) {
191957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
192057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1921f96d37fcSHong Zhang   } else {
1922f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1923f96d37fcSHong Zhang   }
1924f96d37fcSHong Zhang #endif
1925f96d37fcSHong Zhang   PetscFunctionReturn(0);
1926f96d37fcSHong Zhang }
1927