xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 5e68f438519dc746a595544be882e2df42749e8e)
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;
29*5e68f438SHong Zhang   PetscBool      flg;
302499ec78SHong Zhang 
312499ec78SHong Zhang   PetscFunctionBegin;
322499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
33e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
346c4ed002SBarry 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);
35e55be12dSHong Zhang 
360fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
37*5e68f438SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr);
380fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
390fc8cf34SHong Zhang 
40*5e68f438SHong Zhang     if (!flg) { /* set default algorithm based on B->cmap->N */
41*5e68f438SHong Zhang       PetscMPIInt size;
42*5e68f438SHong Zhang       MatInfo     Ainfo,Binfo;
43*5e68f438SHong Zhang       PetscInt    nz_local,nz;
44*5e68f438SHong Zhang 
45*5e68f438SHong Zhang       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
46*5e68f438SHong Zhang       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
47*5e68f438SHong Zhang       ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr);
48*5e68f438SHong Zhang       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
49*5e68f438SHong Zhang       ierr = MPIU_Allreduce(&nz_local,&nz,1,MPIU_INT,MPIU_SUM,comm);CHKERRQ(ierr);
50*5e68f438SHong Zhang 
51*5e68f438SHong Zhang       if (B->cmap->N >100000 && B->cmap->N > fill*(PetscReal)(nz)/size) alg = 0; /* scalable algorithm would take 2x of time as nonscalable algorithm */
52*5e68f438SHong Zhang       ierr = PetscInfo3(B,"Use algorithm %D, BN %D, fill*nz_allocated(A + B)/size %g\n",alg,B->cmap->N,fill*(PetscReal)(nz)/size);CHKERRQ(ierr);
53*5e68f438SHong Zhang     }
54*5e68f438SHong Zhang 
553ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
560fc8cf34SHong Zhang     switch (alg) {
570fc8cf34SHong Zhang     case 1:
580fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
590fc8cf34SHong Zhang       break;
605e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
615e5acdf2Sstefano_zampini     case 2:
625e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
635e5acdf2Sstefano_zampini       break;
645e5acdf2Sstefano_zampini #endif
650fc8cf34SHong Zhang     default:
66b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
670fc8cf34SHong Zhang       break;
680fc8cf34SHong Zhang     }
693ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
702499ec78SHong Zhang   }
713ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
72598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
733ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
742499ec78SHong Zhang   PetscFunctionReturn(0);
752499ec78SHong Zhang }
762499ec78SHong Zhang 
77a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
78598bc09dSHong Zhang {
79598bc09dSHong Zhang   PetscErrorCode ierr;
80598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
81598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
82598bc09dSHong Zhang 
83598bc09dSHong Zhang   PetscFunctionBegin;
84b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
85598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
86a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
87a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
88ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
89a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
90a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
91598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
92598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
93598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
94598bc09dSHong Zhang   PetscFunctionReturn(0);
95598bc09dSHong Zhang }
96598bc09dSHong Zhang 
97a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
984ae0847bSHong Zhang {
994ae0847bSHong Zhang   PetscErrorCode ierr;
1004ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
1014ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
1024ae0847bSHong Zhang 
1034ae0847bSHong Zhang   PetscFunctionBegin;
1044ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
10526fbe8dcSKarl Rupp 
1064ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
1074ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
1084ae0847bSHong Zhang   PetscFunctionReturn(0);
1094ae0847bSHong Zhang }
1104ae0847bSHong Zhang 
111f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
112598bc09dSHong Zhang {
113598bc09dSHong Zhang   PetscErrorCode ierr;
1144ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
115598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
116598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
117c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
118598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
1199ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
120c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
121598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
122c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
12329825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
124598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1259467f45fSHong Zhang   MPI_Comm       comm;
1269467f45fSHong Zhang   PetscMPIInt    size;
127598bc09dSHong Zhang 
128598bc09dSHong Zhang   PetscFunctionBegin;
1299467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1309467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1319467f45fSHong Zhang 
132a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
133598bc09dSHong Zhang   /*-----------------------------------------------------*/
134a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
135b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
136a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
137598bc09dSHong Zhang 
138598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
139598bc09dSHong Zhang   /*----------------------------------------------------------*/
140598bc09dSHong Zhang   /* get data from symbolic products */
141a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
142c12557a1SHong Zhang   p_oth = NULL;
1439467f45fSHong Zhang   if (size >1) {
1449467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1459467f45fSHong Zhang   }
146598bc09dSHong Zhang 
147598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
148598bc09dSHong Zhang   apa = ptap->apa;
149598bc09dSHong Zhang 
15029825010SHong Zhang   api = ptap->api;
15129825010SHong Zhang   apj = ptap->apj;
152598bc09dSHong Zhang   for (i=0; i<cm; i++) {
153c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
154e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
155598bc09dSHong Zhang 
156598bc09dSHong Zhang     /* set values in C */
157598bc09dSHong Zhang     apJ  = apj + api[i];
158598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
159598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
160598bc09dSHong Zhang 
161598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
162598bc09dSHong Zhang     ca = coa + co->i[i];
163598bc09dSHong Zhang     k  = 0;
164598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
165598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
166598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
1675cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
168598bc09dSHong Zhang     }
169598bc09dSHong Zhang 
170598bc09dSHong Zhang     /* diagonal part of C */
171598bc09dSHong Zhang     ca = cda + cd->i[i];
172598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
173598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
1745cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
175598bc09dSHong Zhang     }
176598bc09dSHong Zhang 
177598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
178598bc09dSHong Zhang     ca = coa + co->i[i];
179598bc09dSHong Zhang     for (; k0<conz; k0++) {
180598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
1815cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
182598bc09dSHong Zhang     }
183598bc09dSHong Zhang   }
184598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
185598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
186598bc09dSHong Zhang   PetscFunctionReturn(0);
187598bc09dSHong Zhang }
188598bc09dSHong Zhang 
1890fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
190598bc09dSHong Zhang {
191598bc09dSHong Zhang   PetscErrorCode     ierr;
192ce94432eSBarry Smith   MPI_Comm           comm;
1939467f45fSHong Zhang   PetscMPIInt        size;
194bfcd1a73SHong Zhang   Mat                Cmpi;
195598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1960298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1974ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
198bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1994ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2004ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
201d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
2025cab4f04SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
203598bc09dSHong Zhang   PetscBT            lnkbt;
204598bc09dSHong Zhang   PetscScalar        *apa;
205bfcd1a73SHong Zhang   PetscReal          afill;
206598bc09dSHong Zhang 
207598bc09dSHong Zhang   PetscFunctionBegin;
208ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2099467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2109467f45fSHong Zhang 
211598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
212b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
213598bc09dSHong Zhang 
214598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
215b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
21619f0e57aSHong Zhang 
217598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
218a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
219598bc09dSHong Zhang 
220a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
221598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2229467f45fSHong Zhang   if (size > 1) {
2239467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
224598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
22520e3dc75SHong Zhang   } else {
22620e3dc75SHong Zhang     p_oth = NULL;
22720e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2289467f45fSHong Zhang   }
229598bc09dSHong Zhang 
230598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
231598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
232854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
233a1a4e84aSHong Zhang   ptap->api = api;
234598bc09dSHong Zhang   api[0]    = 0;
235598bc09dSHong Zhang 
2365cab4f04SHong Zhang   /* create and initialize a linked list */
2375cab4f04SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
238598bc09dSHong Zhang 
239bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
240f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
241598bc09dSHong Zhang   current_space = free_space;
242598bc09dSHong Zhang 
243598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
244598bc09dSHong Zhang   for (i=0; i<am; i++) {
245598bc09dSHong Zhang     /* diagonal portion of A */
246598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
247598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
248598bc09dSHong Zhang       row  = *adj++;
249598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
250598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
251598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2521fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
253598bc09dSHong Zhang     }
254598bc09dSHong Zhang     /* off-diagonal portion of A */
255598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
256598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
257598bc09dSHong Zhang       row  = *aoj++;
258598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
259598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2601fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
261598bc09dSHong Zhang     }
262598bc09dSHong Zhang 
263d14fa243SHong Zhang     apnz     = lnk[0];
264598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
265598bc09dSHong Zhang 
266598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
267598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
268f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
269598bc09dSHong Zhang       nspacedouble++;
270598bc09dSHong Zhang     }
271598bc09dSHong Zhang 
272598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2731fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
274598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2752205254eSKarl Rupp 
276598bc09dSHong Zhang     current_space->array           += apnz;
277598bc09dSHong Zhang     current_space->local_used      += apnz;
278598bc09dSHong Zhang     current_space->local_remaining -= apnz;
279598bc09dSHong Zhang   }
280598bc09dSHong Zhang 
281598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
282598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
283854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
284a1a4e84aSHong Zhang   apj  = ptap->apj;
285a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
286598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
287598bc09dSHong Zhang 
288f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2891795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2902205254eSKarl Rupp 
291f84c536bSHong Zhang   ptap->apa = apa;
292f84c536bSHong Zhang 
293bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
294598bc09dSHong Zhang   /*----------------------------------------------------*/
295bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
296bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
29733d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
298a2f3521dSMark F. Adams 
299bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
300bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
301598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
302598bc09dSHong Zhang   for (i=0; i<am; i++) {
303598bc09dSHong Zhang     row  = i + rstart;
304598bc09dSHong Zhang     apnz = api[i+1] - api[i];
305bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
306598bc09dSHong Zhang     apj += apnz;
307598bc09dSHong Zhang   }
308bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
309bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
310598bc09dSHong Zhang 
311bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
312bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
313f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
314bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
315bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
316598bc09dSHong Zhang 
317bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
318bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
319598bc09dSHong Zhang   c->ptap = ptap;
320598bc09dSHong Zhang 
321bfcd1a73SHong Zhang   *C = Cmpi;
322bfcd1a73SHong Zhang 
323bfcd1a73SHong Zhang   /* set MatInfo */
324118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
325bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
326bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
327bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
328bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
329bfcd1a73SHong Zhang 
330bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
331bfcd1a73SHong Zhang   if (api[am]) {
33257622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
33357622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
334bfcd1a73SHong Zhang   } else {
335bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
336bfcd1a73SHong Zhang   }
337bfcd1a73SHong Zhang #endif
338598bc09dSHong Zhang   PetscFunctionReturn(0);
339598bc09dSHong Zhang }
340598bc09dSHong Zhang 
341150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3429123193aSHong Zhang {
3439123193aSHong Zhang   PetscErrorCode ierr;
3449123193aSHong Zhang 
3459123193aSHong Zhang   PetscFunctionBegin;
3469123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3473ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3489123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3493ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3509123193aSHong Zhang   }
3513ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3529123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3533ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3549123193aSHong Zhang   PetscFunctionReturn(0);
3559123193aSHong Zhang }
3569123193aSHong Zhang 
3574324174eSBarry Smith typedef struct {
3584324174eSBarry Smith   Mat         workB;
3594324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3604324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3614324174eSBarry Smith } MPIAIJ_MPIDense;
3624324174eSBarry Smith 
36396b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3644324174eSBarry Smith {
3654324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3664324174eSBarry Smith   PetscErrorCode  ierr;
3674324174eSBarry Smith 
3684324174eSBarry Smith   PetscFunctionBegin;
3696bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3704324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3714324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3724324174eSBarry Smith   PetscFunctionReturn(0);
3734324174eSBarry Smith }
3744324174eSBarry Smith 
375fd4e9aacSBarry Smith /*
376fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
377fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
378fd4e9aacSBarry Smith 
379fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
380fd4e9aacSBarry Smith */
381fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
382fd4e9aacSBarry Smith {
383fd4e9aacSBarry Smith   PetscErrorCode         ierr;
384fd4e9aacSBarry Smith   PetscBool              flg;
385fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
386fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
387fd4e9aacSBarry Smith   PetscContainer         container;
388fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
389fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
390fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
391fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
392fd4e9aacSBarry Smith 
393fd4e9aacSBarry Smith   PetscFunctionBegin;
394fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
395fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
396fd4e9aacSBarry Smith 
397fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
398fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
399fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
400fd4e9aacSBarry Smith 
401fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
402fd4e9aacSBarry Smith 
403fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
404fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
405fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
406fd4e9aacSBarry Smith   /* Create work arrays needed */
407fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
408fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
409fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
410fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
411fd4e9aacSBarry Smith 
412fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
413fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
414fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
415fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
416fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
417fd4e9aacSBarry Smith 
418fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
419fd4e9aacSBarry Smith   PetscFunctionReturn(0);
420fd4e9aacSBarry Smith }
421fd4e9aacSBarry Smith 
4229123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4239123193aSHong Zhang {
4249123193aSHong Zhang   PetscErrorCode         ierr;
425f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
426d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
427bf0cc555SLisandro Dalcin   PetscContainer         container;
4284324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4294324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4304324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4314324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
432d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4339123193aSHong Zhang 
4349123193aSHong Zhang   PetscFunctionBegin;
435ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
436d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
43733d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
438cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4390298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
440cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
441cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4422205254eSKarl Rupp 
443d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
444f32d5d43SBarry Smith 
445b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
446f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4470298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4484324174eSBarry Smith   /* Create work arrays needed */
449dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
450dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
451dcca6d9dSJed Brown                       from->n,&contents->rwaits,
452dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4534324174eSBarry Smith 
454ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
455bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
45696b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
457bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
458bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4599123193aSHong Zhang   PetscFunctionReturn(0);
4609123193aSHong Zhang }
4619123193aSHong Zhang 
462f32d5d43SBarry Smith /*
4632636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4642636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
465f32d5d43SBarry Smith */
4664324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
467f32d5d43SBarry Smith {
468f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
469f32d5d43SBarry Smith   PetscErrorCode         ierr;
470f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
471f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
472f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
473f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
474f32d5d43SBarry Smith   PetscInt               i,j,k;
475aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
476aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
477f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
478ce94432eSBarry Smith   MPI_Comm               comm;
479d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
480f32d5d43SBarry Smith   MPI_Status             status;
4814324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
482bf0cc555SLisandro Dalcin   PetscContainer         container;
4834324174eSBarry Smith   Mat                    workB;
484f32d5d43SBarry Smith 
485f32d5d43SBarry Smith   PetscFunctionBegin;
486ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
487bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
48829825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
489bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4904324174eSBarry Smith 
4914324174eSBarry Smith   workB = *outworkB = contents->workB;
492cf1a0672SHong 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);
493f32d5d43SBarry Smith   sindices = to->indices;
494f32d5d43SBarry Smith   sstarts  = to->starts;
495f32d5d43SBarry Smith   sprocs   = to->procs;
4964324174eSBarry Smith   swaits   = contents->swaits;
4974324174eSBarry Smith   svalues  = contents->svalues;
498f32d5d43SBarry Smith 
499f32d5d43SBarry Smith   rindices = from->indices;
500f32d5d43SBarry Smith   rstarts  = from->starts;
501f32d5d43SBarry Smith   rprocs   = from->procs;
5024324174eSBarry Smith   rwaits   = contents->rwaits;
5034324174eSBarry Smith   rvalues  = contents->rvalues;
504f32d5d43SBarry Smith 
5058c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5068c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
507f32d5d43SBarry Smith 
508f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
509f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
510f32d5d43SBarry Smith   }
511f32d5d43SBarry Smith 
512f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
513f32d5d43SBarry Smith     /* pack a message at a time */
514f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
515f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
516f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5172636ff26SBarry Smith       }
5182636ff26SBarry Smith     }
519f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
520f32d5d43SBarry Smith   }
5212636ff26SBarry Smith 
522f32d5d43SBarry Smith   nrecvs = from->n;
523f32d5d43SBarry Smith   while (nrecvs) {
524f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
525f32d5d43SBarry Smith     nrecvs--;
526f32d5d43SBarry Smith     /* unpack a message at a time */
527f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
528f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
529f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5302636ff26SBarry Smith       }
5312636ff26SBarry Smith     }
532f32d5d43SBarry Smith   }
533cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
534f32d5d43SBarry Smith 
5358c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5368c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
537f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
538f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
539f32d5d43SBarry Smith   PetscFunctionReturn(0);
540f32d5d43SBarry Smith }
541f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
542f32d5d43SBarry Smith 
5439123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5449123193aSHong Zhang {
5459123193aSHong Zhang   PetscErrorCode ierr;
546f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
547f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
548f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
549f32d5d43SBarry Smith   Mat            workB;
5509123193aSHong Zhang 
5519123193aSHong Zhang   PetscFunctionBegin;
552f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
553f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
554f32d5d43SBarry Smith 
555f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5564324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
557f32d5d43SBarry Smith 
558f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
559f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5609123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5619123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5629123193aSHong Zhang   PetscFunctionReturn(0);
5639123193aSHong Zhang }
564cf1a0672SHong Zhang 
565f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5661634b032SHong Zhang {
5671634b032SHong Zhang   PetscErrorCode ierr;
568cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
569cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
570cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
571cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
572cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
573cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
574cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
57529825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
57629825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
577cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
57829825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
579cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
58029825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
58129825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
582a7c7454dSHong Zhang   MPI_Comm       comm;
583a7c7454dSHong Zhang   PetscMPIInt    size;
5841634b032SHong Zhang 
5851634b032SHong Zhang   PetscFunctionBegin;
586a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
587a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
588a7c7454dSHong Zhang 
589cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
590cf1a0672SHong Zhang   /*-----------------------------------------------------*/
591cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
592b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
593cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
594cf1a0672SHong Zhang 
595cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
596cf1a0672SHong Zhang   /*----------------------------------------------------------*/
597cf1a0672SHong Zhang   /* get data from symbolic products */
598cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
599cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
600a7c7454dSHong Zhang   if (size >1) {
601a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
602cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
60320e3dc75SHong Zhang   } else {
60420e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
605a7c7454dSHong Zhang   }
606cf1a0672SHong Zhang 
60729825010SHong Zhang   api = ptap->api;
60829825010SHong Zhang   apj = ptap->apj;
609cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
61029825010SHong Zhang     apJ = apj + api[i];
61129825010SHong Zhang 
612cf1a0672SHong Zhang     /* diagonal portion of A */
613cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
614cf1a0672SHong Zhang     adj = ad->j + adi[i];
615cf1a0672SHong Zhang     ada = ad->a + adi[i];
616cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
617cf1a0672SHong Zhang       row = adj[j];
618cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
619cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
620cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
62129825010SHong Zhang       /* perform sparse axpy */
622cf1a0672SHong Zhang       valtmp = ada[j];
62329825010SHong Zhang       nextp  = 0;
62429825010SHong Zhang       for (k=0; nextp<pnz; k++) {
62529825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
62629825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
62729825010SHong Zhang         }
6281634b032SHong Zhang       }
629cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
630cf1a0672SHong Zhang     }
6311634b032SHong Zhang 
632cf1a0672SHong Zhang     /* off-diagonal portion of A */
633cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
634cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
635cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
636cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
637cf1a0672SHong Zhang       row = aoj[j];
638cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
639cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
640cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
64129825010SHong Zhang       /* perform sparse axpy */
642cf1a0672SHong Zhang       valtmp = aoa[j];
64329825010SHong Zhang       nextp  = 0;
64429825010SHong Zhang       for (k=0; nextp<pnz; k++) {
64529825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
64629825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
64729825010SHong Zhang         }
648cf1a0672SHong Zhang       }
649cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
650cf1a0672SHong Zhang     }
6511634b032SHong Zhang 
652cf1a0672SHong Zhang     /* set values in C */
653cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
654cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6551634b032SHong Zhang 
656cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
657cf1a0672SHong Zhang     ca = coa + co->i[i];
658cf1a0672SHong Zhang     k  = 0;
659cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
660cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
66129825010SHong Zhang       ca[k0]        = apa_sparse[k];
66229825010SHong Zhang       apa_sparse[k] = 0.0;
663cf1a0672SHong Zhang       k++;
664cf1a0672SHong Zhang     }
6651634b032SHong Zhang 
666cf1a0672SHong Zhang     /* diagonal part of C */
667cf1a0672SHong Zhang     ca = cda + cd->i[i];
668cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
66929825010SHong Zhang       ca[k1]        = apa_sparse[k];
67029825010SHong Zhang       apa_sparse[k] = 0.0;
671cf1a0672SHong Zhang       k++;
672cf1a0672SHong Zhang     }
673cf1a0672SHong Zhang 
674cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
675cf1a0672SHong Zhang     ca = coa + co->i[i];
676cf1a0672SHong Zhang     for (; k0<conz; k0++) {
67729825010SHong Zhang       ca[k0]        = apa_sparse[k];
67829825010SHong Zhang       apa_sparse[k] = 0.0;
679cf1a0672SHong Zhang       k++;
680cf1a0672SHong Zhang     }
681cf1a0672SHong Zhang   }
682cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
683cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
684cf1a0672SHong Zhang   PetscFunctionReturn(0);
685cf1a0672SHong Zhang }
686cf1a0672SHong Zhang 
6870fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
688b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
689cf1a0672SHong Zhang {
690cf1a0672SHong Zhang   PetscErrorCode     ierr;
691ce94432eSBarry Smith   MPI_Comm           comm;
692a7c7454dSHong Zhang   PetscMPIInt        size;
693cf1a0672SHong Zhang   Mat                Cmpi;
694cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6950298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
696cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
697cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
698cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
699cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
7000ca7d551SHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max;
7010ca7d551SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
702cf1a0672SHong Zhang   PetscReal          afill;
703cf1a0672SHong Zhang   PetscScalar        *apa;
704eca6b86aSHong Zhang   PetscTable         ta;
705cf1a0672SHong Zhang 
706cf1a0672SHong Zhang   PetscFunctionBegin;
707ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
708a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
709a7c7454dSHong Zhang 
710cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
711b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
712cf1a0672SHong Zhang 
713cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
714b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
71519f0e57aSHong Zhang 
716cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
717cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
718cf1a0672SHong Zhang 
719cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
720cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
721a7c7454dSHong Zhang   if (size > 1) {
722a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
72320e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
724968382fdSHong Zhang   } else {
72520e3dc75SHong Zhang     p_oth  = NULL;
72620e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
727a7c7454dSHong Zhang   }
728cf1a0672SHong Zhang 
729cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
730cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
731854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
732cf1a0672SHong Zhang   ptap->api = api;
733cf1a0672SHong Zhang   api[0]    = 0;
734cf1a0672SHong Zhang 
735cf1a0672SHong Zhang   /* create and initialize a linked list */
73645d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr);
7370ca7d551SHong Zhang 
7380ca7d551SHong Zhang   /* Calculate apnz_max */
7390ca7d551SHong Zhang   apnz_max = 0;
7400ca7d551SHong Zhang   for (i=0; i<am; i++) {
7410ca7d551SHong Zhang     ierr = PetscTableRemoveAll(ta);CHKERRQ(ierr);
7420ca7d551SHong Zhang     /* diagonal portion of A */
7430ca7d551SHong Zhang     nzi  = adi[i+1] - adi[i];
7440ca7d551SHong Zhang     Jptr = adj+adi[i];  /* cols of A_diag */
7450ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_loc,nzi,Jptr,ta);
7460ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7470ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7480ca7d551SHong Zhang 
7490ca7d551SHong Zhang     /*  off-diagonal portion of A */
7500ca7d551SHong Zhang     nzi = aoi[i+1] - aoi[i];
7510ca7d551SHong Zhang     Jptr = aoj+aoi[i];  /* cols of A_off */
7520ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_oth,nzi,Jptr,ta);
7530ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7540ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7550ca7d551SHong Zhang   }
756eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
757eca6b86aSHong Zhang 
7580ca7d551SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(apnz_max,&lnk);CHKERRQ(ierr);
759cf1a0672SHong Zhang 
760cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
761f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
762cf1a0672SHong Zhang   current_space = free_space;
763cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
764cf1a0672SHong Zhang   for (i=0; i<am; i++) {
765cf1a0672SHong Zhang     /* diagonal portion of A */
766cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
767cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
768cf1a0672SHong Zhang       row  = *adj++;
769cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
770cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
771cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
772f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
773cf1a0672SHong Zhang     }
774cf1a0672SHong Zhang     /* off-diagonal portion of A */
775cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
776cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
777cf1a0672SHong Zhang       row  = *aoj++;
778cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
779cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
780f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
781cf1a0672SHong Zhang     }
782cf1a0672SHong Zhang 
783f84c536bSHong Zhang     apnz     = *lnk;
784cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
785cf1a0672SHong Zhang 
786cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
787cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
788f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
789cf1a0672SHong Zhang       nspacedouble++;
790cf1a0672SHong Zhang     }
791cf1a0672SHong Zhang 
792cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
793f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
794cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7952205254eSKarl Rupp 
796cf1a0672SHong Zhang     current_space->array           += apnz;
797cf1a0672SHong Zhang     current_space->local_used      += apnz;
798cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
799cf1a0672SHong Zhang   }
800cf1a0672SHong Zhang 
801cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
802cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
803854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
804cf1a0672SHong Zhang   apj  = ptap->apj;
805cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
806f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
807cf1a0672SHong Zhang 
808cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
809cf1a0672SHong Zhang   /*----------------------------------------------------*/
810cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
811cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
81233d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
813cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
814cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
815cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
816cf1a0672SHong Zhang 
81729825010SHong Zhang   /* malloc apa for assembly Cmpi */
8181795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8192205254eSKarl Rupp 
82029825010SHong Zhang   ptap->apa = apa;
821cf1a0672SHong Zhang   for (i=0; i<am; i++) {
822cf1a0672SHong Zhang     row  = i + rstart;
823cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
824cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
825cf1a0672SHong Zhang     apj += apnz;
826cf1a0672SHong Zhang   }
827cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
828cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
829cf1a0672SHong Zhang 
830cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
831cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
832f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
833cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
834cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
835cf1a0672SHong Zhang 
836cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
837cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
838cf1a0672SHong Zhang   c->ptap = ptap;
839cf1a0672SHong Zhang 
840cf1a0672SHong Zhang   *C = Cmpi;
841cf1a0672SHong Zhang 
842cf1a0672SHong Zhang   /* set MatInfo */
843118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
844cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
845cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
846cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
847cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
848cf1a0672SHong Zhang 
849cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
850cf1a0672SHong Zhang   if (api[am]) {
85157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
85257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
853cf1a0672SHong Zhang   } else {
854cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
855cf1a0672SHong Zhang   }
856cf1a0672SHong Zhang #endif
8571634b032SHong Zhang   PetscFunctionReturn(0);
8581634b032SHong Zhang }
859f96d37fcSHong Zhang 
860f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
861c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
862f96d37fcSHong Zhang {
863f96d37fcSHong Zhang   PetscErrorCode ierr;
864c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
865c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
866f96d37fcSHong Zhang 
867f96d37fcSHong Zhang   PetscFunctionBegin;
868c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
869c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
870c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
871c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
872c216dbf3SHong Zhang 
873c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
874c216dbf3SHong Zhang     switch (alg) {
875c216dbf3SHong Zhang     case 1:
8762bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
877c216dbf3SHong Zhang       break;
878c216dbf3SHong Zhang     case 2:
879275476c6SMatthew G. Knepley     {
880ab1b013aSHong Zhang       Mat         Pt;
881ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
882ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
883ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
884ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
885ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
886ab1b013aSHong Zhang       ptap     = c->ptap;
887ab1b013aSHong Zhang       ptap->Pt = Pt;
888c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
889c216dbf3SHong Zhang       PetscFunctionReturn(0);
890275476c6SMatthew G. Knepley     }
891c216dbf3SHong Zhang       break;
892c216dbf3SHong Zhang     default:
8936da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
894c216dbf3SHong Zhang       break;
895c216dbf3SHong Zhang     }
896c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
897c216dbf3SHong Zhang   }
898c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
899c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
900c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
901ab1b013aSHong Zhang   PetscFunctionReturn(0);
902ab1b013aSHong Zhang }
903ab1b013aSHong Zhang 
904c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
905c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
906c216dbf3SHong Zhang {
907c216dbf3SHong Zhang   PetscErrorCode ierr;
9082bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9092bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9102bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
911c216dbf3SHong Zhang 
912c216dbf3SHong Zhang   PetscFunctionBegin;
913c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
914c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
915f96d37fcSHong Zhang   PetscFunctionReturn(0);
916f96d37fcSHong Zhang }
917f96d37fcSHong Zhang 
9186da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
9196da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
920f96d37fcSHong Zhang {
921c5af039cSHong Zhang   PetscErrorCode      ierr;
922c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
923497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
924c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
925c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
926d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
927497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
928e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
929c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
930ce94432eSBarry Smith   MPI_Comm            comm;
931c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
932c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
933c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
934c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
935c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
936c5af039cSHong Zhang   MPI_Status          *status;
937c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
938d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
939a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
940c5af039cSHong Zhang   Mat                 A_loc;
941c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
942f96d37fcSHong Zhang 
943f96d37fcSHong Zhang   PetscFunctionBegin;
944ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
945c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
946c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
947c5af039cSHong Zhang 
948c5af039cSHong Zhang   ptap  = c->ptap;
949c5af039cSHong Zhang   merge = ptap->merge;
950c5af039cSHong Zhang 
951c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
952c5af039cSHong Zhang   /*--------------------------------------------------------------*/
953c5af039cSHong Zhang   /* get data from symbolic products */
954c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
955854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
956c5af039cSHong Zhang 
957c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
958c5af039cSHong Zhang   owners = merge->rowmap->range;
959854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
960c5af039cSHong Zhang 
961c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
962c5af039cSHong Zhang   A_loc = ptap->A_loc;
963c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
964c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
965d6ab1dc0SHong Zhang   ai    = a_loc->i;
966d6ab1dc0SHong Zhang   aj    = a_loc->j;
967c5af039cSHong Zhang 
968854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
969c5af039cSHong Zhang 
970c5af039cSHong Zhang   for (i=0; i<am; i++) {
971e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
972d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
973d6ab1dc0SHong Zhang     adj = aj + ai[i];
974d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
975c5af039cSHong Zhang     for (j=0; j<anz; j++) {
976e745005bSHong Zhang       aval[adj[j]] = ada[j];
977c5af039cSHong Zhang     }
978c5af039cSHong Zhang 
979c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
980c5af039cSHong Zhang     /*--------------------------------------------------------------*/
981c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
982c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
983c5af039cSHong Zhang     poJ = po->j + po->i[i];
984c5af039cSHong Zhang     pA  = po->a + po->i[i];
985c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
986c5af039cSHong Zhang       row = poJ[j];
987c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
988c5af039cSHong Zhang       cj  = coj + coi[row];
989c5af039cSHong Zhang       ca  = coa + coi[row];
990c5af039cSHong Zhang       /* perform dense axpy */
991c5af039cSHong Zhang       valtmp = pA[j];
992c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
993e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
994c5af039cSHong Zhang       }
995c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
996c5af039cSHong Zhang     }
997c5af039cSHong Zhang 
998c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
999c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1000c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1001c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1002c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1003c5af039cSHong Zhang       row = pdJ[j];
1004c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1005c5af039cSHong Zhang       cj  = bj + bi[row];
1006c5af039cSHong Zhang       ca  = ba + bi[row];
1007c5af039cSHong Zhang       /* perform dense axpy */
1008c5af039cSHong Zhang       valtmp = pA[j];
1009c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1010e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1011c5af039cSHong Zhang       }
1012c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1013c5af039cSHong Zhang     }
1014c5af039cSHong Zhang 
1015d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1016d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1017e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1018c5af039cSHong Zhang   }
1019c5af039cSHong Zhang 
1020c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1021c5af039cSHong Zhang   /*------------------------------------*/
1022c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1023c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1024c5af039cSHong Zhang   len_s  = merge->len_s;
1025c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1026c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1027c5af039cSHong Zhang 
1028dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1029c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1030c5af039cSHong Zhang     if (!len_s[proc]) continue;
1031c5af039cSHong Zhang     i    = merge->owners_co[proc];
1032c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1033c5af039cSHong Zhang     k++;
1034c5af039cSHong Zhang   }
1035c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1036c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1037c5af039cSHong Zhang 
1038c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1039c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1040c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1041c5af039cSHong Zhang 
1042c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1043c5af039cSHong Zhang   /*----------------------------------------------------*/
1044dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1045c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1046c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1047c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1048c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1049c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1050c5af039cSHong Zhang   }
1051c5af039cSHong Zhang 
1052c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1053c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1054c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1055c5af039cSHong Zhang     ba_i = ba + bi[i];
1056c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1057c5af039cSHong Zhang     /* add received vals into ba */
1058c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1059c5af039cSHong Zhang       /* i-th row */
1060c5af039cSHong Zhang       if (i == *nextrow[k]) {
1061c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1062c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1063c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1064c5af039cSHong Zhang         nextcj = 0;
1065c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1066c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1067c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1068c5af039cSHong Zhang           }
1069c5af039cSHong Zhang         }
1070c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1071c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1072c5af039cSHong Zhang       }
1073c5af039cSHong Zhang     }
1074c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1075c5af039cSHong Zhang   }
1076c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1077c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1078c5af039cSHong Zhang 
1079c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1080c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1081c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1082c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1083e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1084f96d37fcSHong Zhang   PetscFunctionReturn(0);
1085f96d37fcSHong Zhang }
1086f96d37fcSHong Zhang 
1087c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1088c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
10892bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1090f96d37fcSHong Zhang {
1091f96d37fcSHong Zhang   PetscErrorCode      ierr;
10924e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1093f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10940298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1095f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1096f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1097f96d37fcSHong Zhang   PetscInt            nnz;
10984e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1099497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1100f96d37fcSHong Zhang   PetscBT             lnkbt;
1101ce94432eSBarry Smith   MPI_Comm            comm;
1102f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1103f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1104f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1105f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1106f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1107f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1108f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1109f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1110d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1111f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1112438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1113f96d37fcSHong Zhang   PetscScalar         *vals;
11144e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1115f96d37fcSHong Zhang 
1116f96d37fcSHong Zhang   PetscFunctionBegin;
1117ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1118c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
11196c4ed002SBarry 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);
1120c5af039cSHong Zhang 
1121f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1122f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1123f96d37fcSHong Zhang 
1124f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1125b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1126f96d37fcSHong Zhang 
1127f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1128f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11292205254eSKarl Rupp 
1130c5af039cSHong Zhang   ptap->A_loc = A_loc;
11312205254eSKarl Rupp 
11321c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1133d6ab1dc0SHong Zhang   ai    = a_loc->i;
1134d6ab1dc0SHong Zhang   aj    = a_loc->j;
1135f96d37fcSHong Zhang 
1136f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1137f96d37fcSHong Zhang   /*----------------------------------------------------*/
11384e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11394e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11404e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11414e938277SHong Zhang 
11424e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11434e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11444e938277SHong Zhang   poti = pot->i; potj = pot->j;
1145f96d37fcSHong Zhang 
1146f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11472205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1148854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1149f96d37fcSHong Zhang   coi[0] = 0;
1150f96d37fcSHong Zhang 
1151f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1152f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1153a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1154f96d37fcSHong Zhang   current_space = free_space;
115519f0e57aSHong Zhang 
1156f96d37fcSHong Zhang   /* create and initialize a linked list */
1157438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1158f96d37fcSHong Zhang 
1159f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1160f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1161f96d37fcSHong Zhang     ptJ = potj + poti[i];
1162f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1163f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1164d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1165d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1166f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1167d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1168f96d37fcSHong Zhang     }
11694e938277SHong Zhang     nnz = lnk[0];
1170f96d37fcSHong Zhang 
1171f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1172f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1173f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1174f96d37fcSHong Zhang       nspacedouble++;
1175f96d37fcSHong Zhang     }
1176f96d37fcSHong Zhang 
1177f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11784e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11792205254eSKarl Rupp 
1180f96d37fcSHong Zhang     current_space->array           += nnz;
1181f96d37fcSHong Zhang     current_space->local_used      += nnz;
1182f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11832205254eSKarl Rupp 
1184f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1185f96d37fcSHong Zhang   }
1186f96d37fcSHong Zhang 
1187854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1188f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11892205254eSKarl Rupp 
1190118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1191f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1192f96d37fcSHong Zhang 
1193f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1194f96d37fcSHong Zhang   /*----------------------------------------------*/
1195f96d37fcSHong Zhang   /* determine row ownership */
1196b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1197f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11982205254eSKarl Rupp 
1199f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1200f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12012205254eSKarl Rupp 
1202f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1203f96d37fcSHong Zhang   owners = merge->rowmap->range;
1204f96d37fcSHong Zhang 
1205f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12061795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1207785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12082205254eSKarl Rupp 
1209f96d37fcSHong Zhang   len_s        = merge->len_s;
1210f96d37fcSHong Zhang   merge->nsend = 0;
1211f96d37fcSHong Zhang 
1212854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1213f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1214f96d37fcSHong Zhang 
1215f96d37fcSHong Zhang   proc = 0;
1216f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1217f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1218f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1219f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1220f96d37fcSHong Zhang   }
1221f96d37fcSHong Zhang 
1222f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1223f96d37fcSHong Zhang   owners_co[0] = 0;
1224f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1225f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1226f96d37fcSHong Zhang     if (len_si[proc]) {
1227f96d37fcSHong Zhang       merge->nsend++;
1228f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1229f96d37fcSHong Zhang       len         += len_si[proc];
1230f96d37fcSHong Zhang     }
1231f96d37fcSHong Zhang   }
1232f96d37fcSHong Zhang 
1233f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12340298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1235f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1236f96d37fcSHong Zhang 
1237f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1238f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1239f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1240854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1241f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1242f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1243f96d37fcSHong Zhang     i    = owners_co[proc];
1244f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1245f96d37fcSHong Zhang     k++;
1246f96d37fcSHong Zhang   }
1247f96d37fcSHong Zhang 
1248f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1249785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1250f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1251c280ed6aSJed Brown     PetscMPIInt icompleted;
1252c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1253f96d37fcSHong Zhang   }
1254f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1255f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1256f96d37fcSHong Zhang 
1257f96d37fcSHong Zhang   /* send and recv coi */
1258f96d37fcSHong Zhang   /*-------------------*/
1259f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1260f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1261854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1262f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1263f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1264f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1265f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1266f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1267f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1268f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1269f96d37fcSHong Zhang     */
1270f96d37fcSHong Zhang     /*-------------------------------------------*/
1271f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1272f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1273f96d37fcSHong Zhang     buf_si[0]   = nrows;
1274f96d37fcSHong Zhang     buf_si_i[0] = 0;
1275f96d37fcSHong Zhang     nrows       = 0;
1276f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1277f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1278f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1279f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1280f96d37fcSHong Zhang       nrows++;
1281f96d37fcSHong Zhang     }
1282f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1283f96d37fcSHong Zhang     k++;
1284f96d37fcSHong Zhang     buf_si += len_si[proc];
1285f96d37fcSHong Zhang   }
1286f96d37fcSHong Zhang   i = merge->nrecv;
1287f96d37fcSHong Zhang   while (i--) {
1288c280ed6aSJed Brown     PetscMPIInt icompleted;
1289c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1290f96d37fcSHong Zhang   }
1291f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1292f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1293f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1294f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1295f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1296f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1297f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1298f96d37fcSHong Zhang 
1299f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1300f96d37fcSHong Zhang   /*------------------------------------------*/
1301f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1302854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1303f96d37fcSHong Zhang   bi[0] = 0;
1304f96d37fcSHong Zhang 
1305c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1306f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1307a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1308f96d37fcSHong Zhang   current_space = free_space;
1309f96d37fcSHong Zhang 
1310dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1311f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1312f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1313f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1314f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1315f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1316f96d37fcSHong Zhang   }
1317f96d37fcSHong Zhang 
13181c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1319f96d37fcSHong Zhang   rmax = 0;
1320f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1321f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1322f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1323f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1324f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1325f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1326d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1327d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1328f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1329d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1330f96d37fcSHong Zhang     }
13314e938277SHong Zhang 
1332f96d37fcSHong Zhang     /* add received col data into lnk */
1333f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1334f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1335f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1336f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13374e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1338f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1339f96d37fcSHong Zhang       }
1340f96d37fcSHong Zhang     }
13414e938277SHong Zhang     nnz = lnk[0];
1342f96d37fcSHong Zhang 
1343f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1344f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1345f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1346f96d37fcSHong Zhang       nspacedouble++;
1347f96d37fcSHong Zhang     }
1348f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13494e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1350f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13512205254eSKarl Rupp 
1352f96d37fcSHong Zhang     current_space->array           += nnz;
1353f96d37fcSHong Zhang     current_space->local_used      += nnz;
1354f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13552205254eSKarl Rupp 
1356f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1357f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1358f96d37fcSHong Zhang   }
1359f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1360f96d37fcSHong Zhang 
1361854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1362f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13632205254eSKarl Rupp 
1364118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1365f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1366d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13674e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13684e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1369f96d37fcSHong Zhang 
13701c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13711c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13721795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1373f96d37fcSHong Zhang 
1374f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1375f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
137633d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1377f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1378f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1379f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1380f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1381f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1382f96d37fcSHong Zhang     row  = i + rstart;
1383f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1384f96d37fcSHong Zhang     Jptr = bj + bi[i];
1385f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1386f96d37fcSHong Zhang   }
1387f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1388f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1389f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1390f96d37fcSHong Zhang 
1391f96d37fcSHong Zhang   merge->bi        = bi;
1392f96d37fcSHong Zhang   merge->bj        = bj;
1393f96d37fcSHong Zhang   merge->coi       = coi;
1394f96d37fcSHong Zhang   merge->coj       = coj;
1395f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1396f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1397f96d37fcSHong Zhang   merge->owners_co = owners_co;
1398f96d37fcSHong Zhang 
1399f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1400f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1401f96d37fcSHong Zhang   c->ptap     = ptap;
14020298fd71SBarry Smith   ptap->api   = NULL;
14030298fd71SBarry Smith   ptap->apj   = NULL;
1404f96d37fcSHong Zhang   ptap->merge = merge;
1405375ed354SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1406375ed354SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1407375ed354SHong Zhang 
1408375ed354SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1409375ed354SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1410375ed354SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1411d6ab1dc0SHong Zhang 
1412d6ab1dc0SHong Zhang   *C = Cmpi;
1413d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1414d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
141557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
141657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1417d6ab1dc0SHong Zhang   } else {
1418d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1419d6ab1dc0SHong Zhang   }
1420d6ab1dc0SHong Zhang #endif
1421d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1422d6ab1dc0SHong Zhang }
1423d6ab1dc0SHong Zhang 
14246da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1425d6ab1dc0SHong Zhang {
1426d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1427d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1428d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1429d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1430d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1431e745005bSHong Zhang   PetscInt            *adj;
1432e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1433e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1434d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1435ce94432eSBarry Smith   MPI_Comm            comm;
1436d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1437d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1438d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1439d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1440d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1441d6ab1dc0SHong Zhang   MPI_Status          *status;
1442d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1443d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1444a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1445d6ab1dc0SHong Zhang   Mat                 A_loc;
1446d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1447d6ab1dc0SHong Zhang 
1448d6ab1dc0SHong Zhang   PetscFunctionBegin;
1449ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1450d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1451d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1452d6ab1dc0SHong Zhang 
1453d6ab1dc0SHong Zhang   ptap  = c->ptap;
1454d6ab1dc0SHong Zhang   merge = ptap->merge;
1455d6ab1dc0SHong Zhang 
1456e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1457e745005bSHong Zhang   /*------------------------------------------*/
1458d6ab1dc0SHong Zhang   /* get data from symbolic products */
1459d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1460854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1461d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1462d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14631795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1464d6ab1dc0SHong Zhang 
1465d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1466d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1467d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1468d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1469d6ab1dc0SHong Zhang   ai    = a_loc->i;
1470d6ab1dc0SHong Zhang   aj    = a_loc->j;
1471d6ab1dc0SHong Zhang 
1472d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1473d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1474d6ab1dc0SHong Zhang     adj = aj + ai[i];
1475d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1476d6ab1dc0SHong Zhang 
1477d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1478e745005bSHong Zhang     /*-------------------------------------------------------------*/
1479d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1480d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1481d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1482d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1483d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1484d6ab1dc0SHong Zhang       row = poJ[j];
1485d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1486d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1487e745005bSHong Zhang       /* perform sparse axpy */
1488e745005bSHong Zhang       nexta  = 0;
1489d6ab1dc0SHong Zhang       valtmp = pA[j];
1490e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1491e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1492e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1493e745005bSHong Zhang           nexta++;
1494d6ab1dc0SHong Zhang         }
1495e745005bSHong Zhang       }
1496e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1497d6ab1dc0SHong Zhang     }
1498d6ab1dc0SHong Zhang 
1499d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1500d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1501d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1502d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1503d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1504d6ab1dc0SHong Zhang       row = pdJ[j];
1505d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1506d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1507e745005bSHong Zhang       /* perform sparse axpy */
1508e745005bSHong Zhang       nexta  = 0;
1509d6ab1dc0SHong Zhang       valtmp = pA[j];
1510e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1511e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1512e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1513e745005bSHong Zhang           nexta++;
1514d6ab1dc0SHong Zhang         }
1515e745005bSHong Zhang       }
1516e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1517d6ab1dc0SHong Zhang     }
1518d6ab1dc0SHong Zhang   }
1519d6ab1dc0SHong Zhang 
1520d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1521d6ab1dc0SHong Zhang   /*------------------------------------*/
1522d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1523d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1524d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1525d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1526d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1527d6ab1dc0SHong Zhang 
1528dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1530d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1531d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1532d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1533d6ab1dc0SHong Zhang     k++;
1534d6ab1dc0SHong Zhang   }
1535d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1536d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1537d6ab1dc0SHong Zhang 
1538d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1539d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1540d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1541d6ab1dc0SHong Zhang 
1542d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1543d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1544dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1545d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1546e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1547d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1548d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1549d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1550d6ab1dc0SHong Zhang   }
1551d6ab1dc0SHong Zhang 
1552d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1553d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1554d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1555d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1556d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1557d6ab1dc0SHong Zhang     /* add received vals into ba */
1558d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1559d6ab1dc0SHong Zhang       /* i-th row */
1560d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1561d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1562d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1563d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1564d6ab1dc0SHong Zhang         nextcj = 0;
1565d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1566d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1567d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1568d6ab1dc0SHong Zhang           }
1569d6ab1dc0SHong Zhang         }
1570d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1571d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1572d6ab1dc0SHong Zhang       }
1573d6ab1dc0SHong Zhang     }
1574d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1575d6ab1dc0SHong Zhang   }
1576d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1577d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1578d6ab1dc0SHong Zhang 
1579d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1580d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1581d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1582d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1583d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1584d6ab1dc0SHong Zhang }
1585d6ab1dc0SHong Zhang 
1586c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
15872bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15882bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
15896da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1590d6ab1dc0SHong Zhang {
1591d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1592d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1593d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15940298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
159533ba5e2fSHong Zhang   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1596d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1597d6ab1dc0SHong Zhang   PetscInt            nnz;
1598d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1599d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1600ce94432eSBarry Smith   MPI_Comm            comm;
1601d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1602d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1603d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1604d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1605d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1606d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1607d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1608d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1609d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1610d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
16114b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1612d6ab1dc0SHong Zhang   PetscScalar         *vals;
1613d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
16140acc65b4SHong Zhang   PetscTable          ta;
1615d6ab1dc0SHong Zhang 
1616d6ab1dc0SHong Zhang   PetscFunctionBegin;
1617ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1618d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
16196c4ed002SBarry 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);
1620d6ab1dc0SHong Zhang 
1621d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1622d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1623d6ab1dc0SHong Zhang 
1624d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1625b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1626d6ab1dc0SHong Zhang 
1627d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1628d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16292205254eSKarl Rupp 
1630d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1631d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1632d6ab1dc0SHong Zhang   ai          = a_loc->i;
1633d6ab1dc0SHong Zhang   aj          = a_loc->j;
1634d6ab1dc0SHong Zhang 
1635d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1636d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1637d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1638d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1639d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1640d6ab1dc0SHong Zhang 
1641d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1642d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1643d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1644d6ab1dc0SHong Zhang 
1645d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1646d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1647d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1648854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1649d6ab1dc0SHong Zhang   coi[0] = 0;
1650d6ab1dc0SHong Zhang 
1651d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1652f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1653a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1654d6ab1dc0SHong Zhang   current_space = free_space;
165519f0e57aSHong Zhang 
1656d6ab1dc0SHong Zhang   /* create and initialize a linked list */
165733ba5e2fSHong Zhang   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
16584b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
16590acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
166033ba5e2fSHong Zhang 
16610acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1662d6ab1dc0SHong Zhang 
1663d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1664d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1665d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1666d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1667d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1668d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1669d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1670d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1671d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1672d6ab1dc0SHong Zhang     }
1673d6ab1dc0SHong Zhang     nnz = lnk[0];
1674d6ab1dc0SHong Zhang 
1675d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1676d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1677f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1678d6ab1dc0SHong Zhang       nspacedouble++;
1679d6ab1dc0SHong Zhang     }
1680d6ab1dc0SHong Zhang 
1681d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1682d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16832205254eSKarl Rupp 
1684d6ab1dc0SHong Zhang     current_space->array           += nnz;
1685d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1686d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16872205254eSKarl Rupp 
1688d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1689d6ab1dc0SHong Zhang   }
1690d6ab1dc0SHong Zhang 
1691854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1692d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16930acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
16942205254eSKarl Rupp 
1695118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1696d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1697d6ab1dc0SHong Zhang 
1698d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1699d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1700d6ab1dc0SHong Zhang   /* determine row ownership */
1701b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1702d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17032205254eSKarl Rupp 
1704d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1705d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17062205254eSKarl Rupp 
1707d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1708d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1709d6ab1dc0SHong Zhang 
1710d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17111795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1712785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17132205254eSKarl Rupp 
1714d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1715d6ab1dc0SHong Zhang   merge->nsend = 0;
1716d6ab1dc0SHong Zhang 
1717854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1718d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1719d6ab1dc0SHong Zhang 
1720d6ab1dc0SHong Zhang   proc = 0;
1721d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1722d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1723d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1724d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1725d6ab1dc0SHong Zhang   }
1726d6ab1dc0SHong Zhang 
1727d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1728d6ab1dc0SHong Zhang   owners_co[0] = 0;
1729d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1730d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1731d6ab1dc0SHong Zhang     if (len_si[proc]) {
1732d6ab1dc0SHong Zhang       merge->nsend++;
1733d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1734d6ab1dc0SHong Zhang       len         += len_si[proc];
1735d6ab1dc0SHong Zhang     }
1736d6ab1dc0SHong Zhang   }
1737d6ab1dc0SHong Zhang 
1738d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17390298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1740d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1741d6ab1dc0SHong Zhang 
1742d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1743d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1745854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1746d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1747d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1748d6ab1dc0SHong Zhang     i    = owners_co[proc];
1749d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1750d6ab1dc0SHong Zhang     k++;
1751d6ab1dc0SHong Zhang   }
1752d6ab1dc0SHong Zhang 
1753d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1754785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1755d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1756c280ed6aSJed Brown     PetscMPIInt icompleted;
1757c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1758d6ab1dc0SHong Zhang   }
1759d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1760d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1761d6ab1dc0SHong Zhang 
17620acc65b4SHong Zhang   /* add received column indices into table to update Armax */
176333ba5e2fSHong Zhang   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
17640acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
17650acc65b4SHong Zhang     Jptr = buf_rj[k];
17660acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
17670acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
17680acc65b4SHong Zhang     }
17690acc65b4SHong Zhang   }
17700acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
177133ba5e2fSHong 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); */
17720acc65b4SHong Zhang 
1773d6ab1dc0SHong Zhang   /* send and recv coi */
1774d6ab1dc0SHong Zhang   /*-------------------*/
1775d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1776d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1777854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1778d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1779d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1780d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1781d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1782d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1783d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1784d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1785d6ab1dc0SHong Zhang     */
1786d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1787d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1788d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1789d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1790d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1791d6ab1dc0SHong Zhang     nrows       = 0;
1792d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1793d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1794d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1795d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1796d6ab1dc0SHong Zhang       nrows++;
1797d6ab1dc0SHong Zhang     }
1798d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1799d6ab1dc0SHong Zhang     k++;
1800d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1801d6ab1dc0SHong Zhang   }
1802d6ab1dc0SHong Zhang   i = merge->nrecv;
1803d6ab1dc0SHong Zhang   while (i--) {
1804c280ed6aSJed Brown     PetscMPIInt icompleted;
1805c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1806d6ab1dc0SHong Zhang   }
1807d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1808d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1809d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1810d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1811d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1812d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1813d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang 
1815d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1816d6ab1dc0SHong Zhang   /*------------------------------------------*/
1817d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1818854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang   bi[0] = 0;
1820d6ab1dc0SHong Zhang 
1821d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1822f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1823a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1824d6ab1dc0SHong Zhang   current_space = free_space;
1825d6ab1dc0SHong Zhang 
1826dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1827d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1828d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1829d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1830d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18312205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1832d6ab1dc0SHong Zhang   }
1833d6ab1dc0SHong Zhang 
18340acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1835d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1836d6ab1dc0SHong Zhang   rmax = 0;
1837d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1838d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1839d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1840d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1841d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1842d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1843d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1844d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1845d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1846d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1847d6ab1dc0SHong Zhang     }
1848d6ab1dc0SHong Zhang 
1849d6ab1dc0SHong Zhang     /* add received col data into lnk */
1850d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1851d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1852d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1853d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1854d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1855d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1856d6ab1dc0SHong Zhang       }
1857d6ab1dc0SHong Zhang     }
1858d6ab1dc0SHong Zhang     nnz = lnk[0];
1859d6ab1dc0SHong Zhang 
1860d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1861d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1862f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1863d6ab1dc0SHong Zhang       nspacedouble++;
1864d6ab1dc0SHong Zhang     }
1865d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1866d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1867d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18682205254eSKarl Rupp 
1869d6ab1dc0SHong Zhang     current_space->array           += nnz;
1870d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1871d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18722205254eSKarl Rupp 
1873d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1874d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1875d6ab1dc0SHong Zhang   }
1876f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1877d6ab1dc0SHong Zhang 
1878854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1879d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1880118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1881d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1882d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
18830acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
18840acc65b4SHong Zhang 
1885d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1886d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1887d6ab1dc0SHong Zhang 
1888d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1889d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18901795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1891d6ab1dc0SHong Zhang 
1892d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1893d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
189433d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1895d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1896d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1897d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1898d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1899d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1900d6ab1dc0SHong Zhang     row  = i + rstart;
1901d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1902d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1903d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1904d6ab1dc0SHong Zhang   }
1905d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1906d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1907d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1908d6ab1dc0SHong Zhang 
1909d6ab1dc0SHong Zhang   merge->bi        = bi;
1910d6ab1dc0SHong Zhang   merge->bj        = bj;
1911d6ab1dc0SHong Zhang   merge->coi       = coi;
1912d6ab1dc0SHong Zhang   merge->coj       = coj;
1913d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1914d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1915d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1916d6ab1dc0SHong Zhang 
1917d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1918d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19192205254eSKarl Rupp 
1920d6ab1dc0SHong Zhang   c->ptap     = ptap;
19210298fd71SBarry Smith   ptap->api   = NULL;
19220298fd71SBarry Smith   ptap->apj   = NULL;
1923d6ab1dc0SHong Zhang   ptap->merge = merge;
19240298fd71SBarry Smith   ptap->apa   = NULL;
1925a560ef98SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1926a560ef98SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1927a560ef98SHong Zhang 
1928a560ef98SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1929a560ef98SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1930a560ef98SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1931f96d37fcSHong Zhang 
1932f96d37fcSHong Zhang   *C = Cmpi;
1933f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1934f96d37fcSHong Zhang   if (bi[pn] != 0) {
193557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
193657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1937f96d37fcSHong Zhang   } else {
1938f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1939f96d37fcSHong Zhang   }
1940f96d37fcSHong Zhang #endif
1941f96d37fcSHong Zhang   PetscFunctionReturn(0);
1942f96d37fcSHong Zhang }
1943