xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 60106aa78d02592cecf2730597e908291462046e)
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
27*60106aa7SHong Zhang   PetscInt       alg = 1; /* set nonscalable algorithm as default */
28e55be12dSHong Zhang   MPI_Comm       comm;
295e68f438SHong 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);
375e68f438SHong 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*60106aa7SHong Zhang     if (!flg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
415e68f438SHong Zhang       MatInfo     Ainfo,Binfo;
42*60106aa7SHong Zhang       PetscInt    nz_local;
43*60106aa7SHong Zhang       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
445e68f438SHong Zhang 
455e68f438SHong Zhang       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
465e68f438SHong Zhang       ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr);
475e68f438SHong Zhang       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
485e68f438SHong Zhang 
49*60106aa7SHong Zhang       if (B->cmap->N > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
50*60106aa7SHong Zhang       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
51*60106aa7SHong Zhang 
52*60106aa7SHong Zhang       if (alg_scalable) {
53*60106aa7SHong Zhang         alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
54*60106aa7SHong Zhang         ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,fill*nz_local);CHKERRQ(ierr);
55*60106aa7SHong Zhang       }
565e68f438SHong Zhang     }
575e68f438SHong Zhang 
583ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
590fc8cf34SHong Zhang     switch (alg) {
600fc8cf34SHong Zhang     case 1:
610fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
620fc8cf34SHong Zhang       break;
635e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
645e5acdf2Sstefano_zampini     case 2:
655e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
665e5acdf2Sstefano_zampini       break;
675e5acdf2Sstefano_zampini #endif
680fc8cf34SHong Zhang     default:
69b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
700fc8cf34SHong Zhang       break;
710fc8cf34SHong Zhang     }
723ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
732499ec78SHong Zhang   }
743ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
75598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
763ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
772499ec78SHong Zhang   PetscFunctionReturn(0);
782499ec78SHong Zhang }
792499ec78SHong Zhang 
80a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
81598bc09dSHong Zhang {
82598bc09dSHong Zhang   PetscErrorCode ierr;
83598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
84598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
85598bc09dSHong Zhang 
86598bc09dSHong Zhang   PetscFunctionBegin;
87b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
88598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
89a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
90a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
91ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
92a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
93a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
94598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
95598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
96598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
97598bc09dSHong Zhang   PetscFunctionReturn(0);
98598bc09dSHong Zhang }
99598bc09dSHong Zhang 
100a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
1014ae0847bSHong Zhang {
1024ae0847bSHong Zhang   PetscErrorCode ierr;
1034ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
1044ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
1054ae0847bSHong Zhang 
1064ae0847bSHong Zhang   PetscFunctionBegin;
1074ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
10826fbe8dcSKarl Rupp 
1094ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
1104ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
1114ae0847bSHong Zhang   PetscFunctionReturn(0);
1124ae0847bSHong Zhang }
1134ae0847bSHong Zhang 
114f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
115598bc09dSHong Zhang {
116598bc09dSHong Zhang   PetscErrorCode ierr;
1174ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
118598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
119598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
120c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
121598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
1229ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
123c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
124598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
125c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
12629825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
127598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1289467f45fSHong Zhang   MPI_Comm       comm;
1299467f45fSHong Zhang   PetscMPIInt    size;
130598bc09dSHong Zhang 
131598bc09dSHong Zhang   PetscFunctionBegin;
1329467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1339467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1349467f45fSHong Zhang 
135a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
136598bc09dSHong Zhang   /*-----------------------------------------------------*/
137a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
138b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
139a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
140598bc09dSHong Zhang 
141598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
142598bc09dSHong Zhang   /*----------------------------------------------------------*/
143598bc09dSHong Zhang   /* get data from symbolic products */
144a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
145c12557a1SHong Zhang   p_oth = NULL;
1469467f45fSHong Zhang   if (size >1) {
1479467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1489467f45fSHong Zhang   }
149598bc09dSHong Zhang 
150598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
151598bc09dSHong Zhang   apa = ptap->apa;
152598bc09dSHong Zhang 
15329825010SHong Zhang   api = ptap->api;
15429825010SHong Zhang   apj = ptap->apj;
155598bc09dSHong Zhang   for (i=0; i<cm; i++) {
156c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
157e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
158598bc09dSHong Zhang 
159598bc09dSHong Zhang     /* set values in C */
160598bc09dSHong Zhang     apJ  = apj + api[i];
161598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
162598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
163598bc09dSHong Zhang 
164598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
165598bc09dSHong Zhang     ca = coa + co->i[i];
166598bc09dSHong Zhang     k  = 0;
167598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
168598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
169598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
1705cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
171598bc09dSHong Zhang     }
172598bc09dSHong Zhang 
173598bc09dSHong Zhang     /* diagonal part of C */
174598bc09dSHong Zhang     ca = cda + cd->i[i];
175598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
176598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
1775cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
178598bc09dSHong Zhang     }
179598bc09dSHong Zhang 
180598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
181598bc09dSHong Zhang     ca = coa + co->i[i];
182598bc09dSHong Zhang     for (; k0<conz; k0++) {
183598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
1845cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
185598bc09dSHong Zhang     }
186598bc09dSHong Zhang   }
187598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
188598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
189598bc09dSHong Zhang   PetscFunctionReturn(0);
190598bc09dSHong Zhang }
191598bc09dSHong Zhang 
1920fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
193598bc09dSHong Zhang {
194598bc09dSHong Zhang   PetscErrorCode     ierr;
195ce94432eSBarry Smith   MPI_Comm           comm;
1969467f45fSHong Zhang   PetscMPIInt        size;
197bfcd1a73SHong Zhang   Mat                Cmpi;
198598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1990298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2004ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
201bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2024ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2034ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
204d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
2055cab4f04SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
206598bc09dSHong Zhang   PetscBT            lnkbt;
207598bc09dSHong Zhang   PetscScalar        *apa;
208bfcd1a73SHong Zhang   PetscReal          afill;
209598bc09dSHong Zhang 
210598bc09dSHong Zhang   PetscFunctionBegin;
211ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2129467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2139467f45fSHong Zhang 
214598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
215b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
216598bc09dSHong Zhang 
217598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
218b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
21919f0e57aSHong Zhang 
220598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
221a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
222598bc09dSHong Zhang 
223a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
224598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2259467f45fSHong Zhang   if (size > 1) {
2269467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
227598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
22820e3dc75SHong Zhang   } else {
22920e3dc75SHong Zhang     p_oth = NULL;
23020e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2319467f45fSHong Zhang   }
232598bc09dSHong Zhang 
233598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
234598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
235854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
236a1a4e84aSHong Zhang   ptap->api = api;
237598bc09dSHong Zhang   api[0]    = 0;
238598bc09dSHong Zhang 
2395cab4f04SHong Zhang   /* create and initialize a linked list */
2405cab4f04SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
241598bc09dSHong Zhang 
242bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
243f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
244598bc09dSHong Zhang   current_space = free_space;
245598bc09dSHong Zhang 
246598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
247598bc09dSHong Zhang   for (i=0; i<am; i++) {
248598bc09dSHong Zhang     /* diagonal portion of A */
249598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
250598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
251598bc09dSHong Zhang       row  = *adj++;
252598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
253598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
254598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2551fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
256598bc09dSHong Zhang     }
257598bc09dSHong Zhang     /* off-diagonal portion of A */
258598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
259598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
260598bc09dSHong Zhang       row  = *aoj++;
261598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
262598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2631fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
264598bc09dSHong Zhang     }
265598bc09dSHong Zhang 
266d14fa243SHong Zhang     apnz     = lnk[0];
267598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
268598bc09dSHong Zhang 
269598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
270598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
271f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
272598bc09dSHong Zhang       nspacedouble++;
273598bc09dSHong Zhang     }
274598bc09dSHong Zhang 
275598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2761fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
277598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2782205254eSKarl Rupp 
279598bc09dSHong Zhang     current_space->array           += apnz;
280598bc09dSHong Zhang     current_space->local_used      += apnz;
281598bc09dSHong Zhang     current_space->local_remaining -= apnz;
282598bc09dSHong Zhang   }
283598bc09dSHong Zhang 
284598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
285598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
286854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
287a1a4e84aSHong Zhang   apj  = ptap->apj;
288a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
289598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
290598bc09dSHong Zhang 
291f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2921795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2932205254eSKarl Rupp 
294f84c536bSHong Zhang   ptap->apa = apa;
295f84c536bSHong Zhang 
296bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
297598bc09dSHong Zhang   /*----------------------------------------------------*/
298bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
299bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
30033d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
301a2f3521dSMark F. Adams 
302bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
303bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
304598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
305598bc09dSHong Zhang   for (i=0; i<am; i++) {
306598bc09dSHong Zhang     row  = i + rstart;
307598bc09dSHong Zhang     apnz = api[i+1] - api[i];
308bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
309598bc09dSHong Zhang     apj += apnz;
310598bc09dSHong Zhang   }
311bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313598bc09dSHong Zhang 
314bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
315bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
316f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
317bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
318bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
319598bc09dSHong Zhang 
320bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
321bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
322598bc09dSHong Zhang   c->ptap = ptap;
323598bc09dSHong Zhang 
324bfcd1a73SHong Zhang   *C = Cmpi;
325bfcd1a73SHong Zhang 
326bfcd1a73SHong Zhang   /* set MatInfo */
327118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
328bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
329bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
330bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
331bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
332bfcd1a73SHong Zhang 
333bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
334bfcd1a73SHong Zhang   if (api[am]) {
33557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
33657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
337bfcd1a73SHong Zhang   } else {
338bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
339bfcd1a73SHong Zhang   }
340bfcd1a73SHong Zhang #endif
341598bc09dSHong Zhang   PetscFunctionReturn(0);
342598bc09dSHong Zhang }
343598bc09dSHong Zhang 
344150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3459123193aSHong Zhang {
3469123193aSHong Zhang   PetscErrorCode ierr;
3479123193aSHong Zhang 
3489123193aSHong Zhang   PetscFunctionBegin;
3499123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3503ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3519123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3523ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3539123193aSHong Zhang   }
3543ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3559123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3563ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3579123193aSHong Zhang   PetscFunctionReturn(0);
3589123193aSHong Zhang }
3599123193aSHong Zhang 
3604324174eSBarry Smith typedef struct {
3614324174eSBarry Smith   Mat         workB;
3624324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3634324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3644324174eSBarry Smith } MPIAIJ_MPIDense;
3654324174eSBarry Smith 
36696b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3674324174eSBarry Smith {
3684324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3694324174eSBarry Smith   PetscErrorCode  ierr;
3704324174eSBarry Smith 
3714324174eSBarry Smith   PetscFunctionBegin;
3726bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3734324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3744324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3754324174eSBarry Smith   PetscFunctionReturn(0);
3764324174eSBarry Smith }
3774324174eSBarry Smith 
378fd4e9aacSBarry Smith /*
379fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
380fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
381fd4e9aacSBarry Smith 
382fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
383fd4e9aacSBarry Smith */
384fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
385fd4e9aacSBarry Smith {
386fd4e9aacSBarry Smith   PetscErrorCode         ierr;
387fd4e9aacSBarry Smith   PetscBool              flg;
388fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
389fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
390fd4e9aacSBarry Smith   PetscContainer         container;
391fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
392fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
393fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
394fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
395fd4e9aacSBarry Smith 
396fd4e9aacSBarry Smith   PetscFunctionBegin;
397fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
398fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
399fd4e9aacSBarry Smith 
400fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
401fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
402fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
403fd4e9aacSBarry Smith 
404fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
405fd4e9aacSBarry Smith 
406fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
407fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
408fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
409fd4e9aacSBarry Smith   /* Create work arrays needed */
410fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
411fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
412fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
413fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
414fd4e9aacSBarry Smith 
415fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
416fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
417fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
418fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
419fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
420fd4e9aacSBarry Smith 
421fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
422fd4e9aacSBarry Smith   PetscFunctionReturn(0);
423fd4e9aacSBarry Smith }
424fd4e9aacSBarry Smith 
4259123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4269123193aSHong Zhang {
4279123193aSHong Zhang   PetscErrorCode         ierr;
428f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
429d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
430bf0cc555SLisandro Dalcin   PetscContainer         container;
4314324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4324324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4334324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4344324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
435d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4369123193aSHong Zhang 
4379123193aSHong Zhang   PetscFunctionBegin;
438ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
439d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
44033d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
441cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4420298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
443cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
444cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4452205254eSKarl Rupp 
446d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
447f32d5d43SBarry Smith 
448b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
449f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4500298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4514324174eSBarry Smith   /* Create work arrays needed */
452dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
453dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
454dcca6d9dSJed Brown                       from->n,&contents->rwaits,
455dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4564324174eSBarry Smith 
457ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
458bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
45996b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
460bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
461bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4629123193aSHong Zhang   PetscFunctionReturn(0);
4639123193aSHong Zhang }
4649123193aSHong Zhang 
465f32d5d43SBarry Smith /*
4662636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4672636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
468f32d5d43SBarry Smith */
4694324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
470f32d5d43SBarry Smith {
471f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
472f32d5d43SBarry Smith   PetscErrorCode         ierr;
473f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
474f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
475f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
476f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
477f32d5d43SBarry Smith   PetscInt               i,j,k;
478aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
479aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
480f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
481ce94432eSBarry Smith   MPI_Comm               comm;
482d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
483f32d5d43SBarry Smith   MPI_Status             status;
4844324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
485bf0cc555SLisandro Dalcin   PetscContainer         container;
4864324174eSBarry Smith   Mat                    workB;
487f32d5d43SBarry Smith 
488f32d5d43SBarry Smith   PetscFunctionBegin;
489ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
490bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
49129825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
492bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4934324174eSBarry Smith 
4944324174eSBarry Smith   workB = *outworkB = contents->workB;
495cf1a0672SHong 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);
496f32d5d43SBarry Smith   sindices = to->indices;
497f32d5d43SBarry Smith   sstarts  = to->starts;
498f32d5d43SBarry Smith   sprocs   = to->procs;
4994324174eSBarry Smith   swaits   = contents->swaits;
5004324174eSBarry Smith   svalues  = contents->svalues;
501f32d5d43SBarry Smith 
502f32d5d43SBarry Smith   rindices = from->indices;
503f32d5d43SBarry Smith   rstarts  = from->starts;
504f32d5d43SBarry Smith   rprocs   = from->procs;
5054324174eSBarry Smith   rwaits   = contents->rwaits;
5064324174eSBarry Smith   rvalues  = contents->rvalues;
507f32d5d43SBarry Smith 
5088c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5098c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
510f32d5d43SBarry Smith 
511f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
512f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
513f32d5d43SBarry Smith   }
514f32d5d43SBarry Smith 
515f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
516f32d5d43SBarry Smith     /* pack a message at a time */
517f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
518f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
519f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5202636ff26SBarry Smith       }
5212636ff26SBarry Smith     }
522f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
523f32d5d43SBarry Smith   }
5242636ff26SBarry Smith 
525f32d5d43SBarry Smith   nrecvs = from->n;
526f32d5d43SBarry Smith   while (nrecvs) {
527f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
528f32d5d43SBarry Smith     nrecvs--;
529f32d5d43SBarry Smith     /* unpack a message at a time */
530f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
531f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
532f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5332636ff26SBarry Smith       }
5342636ff26SBarry Smith     }
535f32d5d43SBarry Smith   }
536cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
537f32d5d43SBarry Smith 
5388c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5398c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
540f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
541f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
542f32d5d43SBarry Smith   PetscFunctionReturn(0);
543f32d5d43SBarry Smith }
544f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
545f32d5d43SBarry Smith 
5469123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5479123193aSHong Zhang {
5489123193aSHong Zhang   PetscErrorCode ierr;
549f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
550f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
551f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
552f32d5d43SBarry Smith   Mat            workB;
5539123193aSHong Zhang 
5549123193aSHong Zhang   PetscFunctionBegin;
555f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
556f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
557f32d5d43SBarry Smith 
558f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5594324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
560f32d5d43SBarry Smith 
561f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
562f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5639123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5649123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5659123193aSHong Zhang   PetscFunctionReturn(0);
5669123193aSHong Zhang }
567cf1a0672SHong Zhang 
568f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5691634b032SHong Zhang {
5701634b032SHong Zhang   PetscErrorCode ierr;
571cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
572cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
573cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
574cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
575cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
576cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
577cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
57829825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
57929825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
580cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
58129825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
582cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
58329825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
58429825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
585a7c7454dSHong Zhang   MPI_Comm       comm;
586a7c7454dSHong Zhang   PetscMPIInt    size;
5871634b032SHong Zhang 
5881634b032SHong Zhang   PetscFunctionBegin;
589a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
590a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
591a7c7454dSHong Zhang 
592cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
593cf1a0672SHong Zhang   /*-----------------------------------------------------*/
594cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
595b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
596cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
597cf1a0672SHong Zhang 
598cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
599cf1a0672SHong Zhang   /*----------------------------------------------------------*/
600cf1a0672SHong Zhang   /* get data from symbolic products */
601cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
602cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
603a7c7454dSHong Zhang   if (size >1) {
604a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
605cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
60620e3dc75SHong Zhang   } else {
60720e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
608a7c7454dSHong Zhang   }
609cf1a0672SHong Zhang 
61029825010SHong Zhang   api = ptap->api;
61129825010SHong Zhang   apj = ptap->apj;
612cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
61329825010SHong Zhang     apJ = apj + api[i];
61429825010SHong Zhang 
615cf1a0672SHong Zhang     /* diagonal portion of A */
616cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
617cf1a0672SHong Zhang     adj = ad->j + adi[i];
618cf1a0672SHong Zhang     ada = ad->a + adi[i];
619cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
620cf1a0672SHong Zhang       row = adj[j];
621cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
622cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
623cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
62429825010SHong Zhang       /* perform sparse axpy */
625cf1a0672SHong Zhang       valtmp = ada[j];
62629825010SHong Zhang       nextp  = 0;
62729825010SHong Zhang       for (k=0; nextp<pnz; k++) {
62829825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
62929825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63029825010SHong Zhang         }
6311634b032SHong Zhang       }
632cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
633cf1a0672SHong Zhang     }
6341634b032SHong Zhang 
635cf1a0672SHong Zhang     /* off-diagonal portion of A */
636cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
637cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
638cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
639cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
640cf1a0672SHong Zhang       row = aoj[j];
641cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
642cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
643cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
64429825010SHong Zhang       /* perform sparse axpy */
645cf1a0672SHong Zhang       valtmp = aoa[j];
64629825010SHong Zhang       nextp  = 0;
64729825010SHong Zhang       for (k=0; nextp<pnz; k++) {
64829825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
64929825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
65029825010SHong Zhang         }
651cf1a0672SHong Zhang       }
652cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
653cf1a0672SHong Zhang     }
6541634b032SHong Zhang 
655cf1a0672SHong Zhang     /* set values in C */
656cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
657cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6581634b032SHong Zhang 
659cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
660cf1a0672SHong Zhang     ca = coa + co->i[i];
661cf1a0672SHong Zhang     k  = 0;
662cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
663cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
66429825010SHong Zhang       ca[k0]        = apa_sparse[k];
66529825010SHong Zhang       apa_sparse[k] = 0.0;
666cf1a0672SHong Zhang       k++;
667cf1a0672SHong Zhang     }
6681634b032SHong Zhang 
669cf1a0672SHong Zhang     /* diagonal part of C */
670cf1a0672SHong Zhang     ca = cda + cd->i[i];
671cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
67229825010SHong Zhang       ca[k1]        = apa_sparse[k];
67329825010SHong Zhang       apa_sparse[k] = 0.0;
674cf1a0672SHong Zhang       k++;
675cf1a0672SHong Zhang     }
676cf1a0672SHong Zhang 
677cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
678cf1a0672SHong Zhang     ca = coa + co->i[i];
679cf1a0672SHong Zhang     for (; k0<conz; k0++) {
68029825010SHong Zhang       ca[k0]        = apa_sparse[k];
68129825010SHong Zhang       apa_sparse[k] = 0.0;
682cf1a0672SHong Zhang       k++;
683cf1a0672SHong Zhang     }
684cf1a0672SHong Zhang   }
685cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
687cf1a0672SHong Zhang   PetscFunctionReturn(0);
688cf1a0672SHong Zhang }
689cf1a0672SHong Zhang 
6900fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
691b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
692cf1a0672SHong Zhang {
693cf1a0672SHong Zhang   PetscErrorCode     ierr;
694ce94432eSBarry Smith   MPI_Comm           comm;
695a7c7454dSHong Zhang   PetscMPIInt        size;
696cf1a0672SHong Zhang   Mat                Cmpi;
697cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6980298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
699cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
700cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
701cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
702cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
7030ca7d551SHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max;
7040ca7d551SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
705cf1a0672SHong Zhang   PetscReal          afill;
706cf1a0672SHong Zhang   PetscScalar        *apa;
707eca6b86aSHong Zhang   PetscTable         ta;
708cf1a0672SHong Zhang 
709cf1a0672SHong Zhang   PetscFunctionBegin;
710ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
711a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
712a7c7454dSHong Zhang 
713cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
714b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
715cf1a0672SHong Zhang 
716cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
717b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
71819f0e57aSHong Zhang 
719cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
720cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
721cf1a0672SHong Zhang 
722cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
723cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
724a7c7454dSHong Zhang   if (size > 1) {
725a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
72620e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
727968382fdSHong Zhang   } else {
72820e3dc75SHong Zhang     p_oth  = NULL;
72920e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
730a7c7454dSHong Zhang   }
731cf1a0672SHong Zhang 
732cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
733cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
734854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
735cf1a0672SHong Zhang   ptap->api = api;
736cf1a0672SHong Zhang   api[0]    = 0;
737cf1a0672SHong Zhang 
738cf1a0672SHong Zhang   /* create and initialize a linked list */
73945d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr);
7400ca7d551SHong Zhang 
7410ca7d551SHong Zhang   /* Calculate apnz_max */
7420ca7d551SHong Zhang   apnz_max = 0;
7430ca7d551SHong Zhang   for (i=0; i<am; i++) {
7440ca7d551SHong Zhang     ierr = PetscTableRemoveAll(ta);CHKERRQ(ierr);
7450ca7d551SHong Zhang     /* diagonal portion of A */
7460ca7d551SHong Zhang     nzi  = adi[i+1] - adi[i];
7470ca7d551SHong Zhang     Jptr = adj+adi[i];  /* cols of A_diag */
7480ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_loc,nzi,Jptr,ta);
7490ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7500ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7510ca7d551SHong Zhang 
7520ca7d551SHong Zhang     /*  off-diagonal portion of A */
7530ca7d551SHong Zhang     nzi = aoi[i+1] - aoi[i];
7540ca7d551SHong Zhang     Jptr = aoj+aoi[i];  /* cols of A_off */
7550ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_oth,nzi,Jptr,ta);
7560ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7570ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7580ca7d551SHong Zhang   }
759eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
760eca6b86aSHong Zhang 
7610ca7d551SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(apnz_max,&lnk);CHKERRQ(ierr);
762cf1a0672SHong Zhang 
763cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
764f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
765cf1a0672SHong Zhang   current_space = free_space;
766cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
767cf1a0672SHong Zhang   for (i=0; i<am; i++) {
768cf1a0672SHong Zhang     /* diagonal portion of A */
769cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
770cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
771cf1a0672SHong Zhang       row  = *adj++;
772cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
773cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
774cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
775f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
776cf1a0672SHong Zhang     }
777cf1a0672SHong Zhang     /* off-diagonal portion of A */
778cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
779cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
780cf1a0672SHong Zhang       row  = *aoj++;
781cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
782cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
783f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
784cf1a0672SHong Zhang     }
785cf1a0672SHong Zhang 
786f84c536bSHong Zhang     apnz     = *lnk;
787cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
788cf1a0672SHong Zhang 
789cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
790cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
791f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
792cf1a0672SHong Zhang       nspacedouble++;
793cf1a0672SHong Zhang     }
794cf1a0672SHong Zhang 
795cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
796f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
797cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7982205254eSKarl Rupp 
799cf1a0672SHong Zhang     current_space->array           += apnz;
800cf1a0672SHong Zhang     current_space->local_used      += apnz;
801cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
802cf1a0672SHong Zhang   }
803cf1a0672SHong Zhang 
804cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
805cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
806854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
807cf1a0672SHong Zhang   apj  = ptap->apj;
808cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
809f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
810cf1a0672SHong Zhang 
811cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
812cf1a0672SHong Zhang   /*----------------------------------------------------*/
813cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
814cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
81533d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
816cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
817cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
818cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
819cf1a0672SHong Zhang 
82029825010SHong Zhang   /* malloc apa for assembly Cmpi */
8211795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8222205254eSKarl Rupp 
82329825010SHong Zhang   ptap->apa = apa;
824cf1a0672SHong Zhang   for (i=0; i<am; i++) {
825cf1a0672SHong Zhang     row  = i + rstart;
826cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
827cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
828cf1a0672SHong Zhang     apj += apnz;
829cf1a0672SHong Zhang   }
830cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
831cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
832cf1a0672SHong Zhang 
833cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
834cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
835f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
836cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
837cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
838cf1a0672SHong Zhang 
839cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
840cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
841cf1a0672SHong Zhang   c->ptap = ptap;
842cf1a0672SHong Zhang 
843cf1a0672SHong Zhang   *C = Cmpi;
844cf1a0672SHong Zhang 
845cf1a0672SHong Zhang   /* set MatInfo */
846118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
847cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
848cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
849cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
850cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
851cf1a0672SHong Zhang 
852cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
853cf1a0672SHong Zhang   if (api[am]) {
85457622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
85557622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
856cf1a0672SHong Zhang   } else {
857cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
858cf1a0672SHong Zhang   }
859cf1a0672SHong Zhang #endif
8601634b032SHong Zhang   PetscFunctionReturn(0);
8611634b032SHong Zhang }
862f96d37fcSHong Zhang 
863f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
864c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
865f96d37fcSHong Zhang {
866f96d37fcSHong Zhang   PetscErrorCode ierr;
867c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
868c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
869f96d37fcSHong Zhang 
870f96d37fcSHong Zhang   PetscFunctionBegin;
871c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
872c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
873c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
874c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
875c216dbf3SHong Zhang 
876c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
877c216dbf3SHong Zhang     switch (alg) {
878c216dbf3SHong Zhang     case 1:
8792bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
880c216dbf3SHong Zhang       break;
881c216dbf3SHong Zhang     case 2:
882275476c6SMatthew G. Knepley     {
883ab1b013aSHong Zhang       Mat         Pt;
884ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
885ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
886ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
887ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
888ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
889ab1b013aSHong Zhang       ptap     = c->ptap;
890ab1b013aSHong Zhang       ptap->Pt = Pt;
891c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
892c216dbf3SHong Zhang       PetscFunctionReturn(0);
893275476c6SMatthew G. Knepley     }
894c216dbf3SHong Zhang       break;
895c216dbf3SHong Zhang     default:
8966da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
897c216dbf3SHong Zhang       break;
898c216dbf3SHong Zhang     }
899c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
900c216dbf3SHong Zhang   }
901c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
902c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
903c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
904ab1b013aSHong Zhang   PetscFunctionReturn(0);
905ab1b013aSHong Zhang }
906ab1b013aSHong Zhang 
907c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
908c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
909c216dbf3SHong Zhang {
910c216dbf3SHong Zhang   PetscErrorCode ierr;
9112bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9122bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9132bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
914c216dbf3SHong Zhang 
915c216dbf3SHong Zhang   PetscFunctionBegin;
916c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
917c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
918f96d37fcSHong Zhang   PetscFunctionReturn(0);
919f96d37fcSHong Zhang }
920f96d37fcSHong Zhang 
9216da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
9226da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
923f96d37fcSHong Zhang {
924c5af039cSHong Zhang   PetscErrorCode      ierr;
925c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
926497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
927c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
928c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
929d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
930497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
931e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
932c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
933ce94432eSBarry Smith   MPI_Comm            comm;
934c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
935c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
936c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
937c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
938c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
939c5af039cSHong Zhang   MPI_Status          *status;
940c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
941d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
942a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
943c5af039cSHong Zhang   Mat                 A_loc;
944c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
945f96d37fcSHong Zhang 
946f96d37fcSHong Zhang   PetscFunctionBegin;
947ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
948c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
949c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
950c5af039cSHong Zhang 
951c5af039cSHong Zhang   ptap  = c->ptap;
952c5af039cSHong Zhang   merge = ptap->merge;
953c5af039cSHong Zhang 
954c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
955c5af039cSHong Zhang   /*--------------------------------------------------------------*/
956c5af039cSHong Zhang   /* get data from symbolic products */
957c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
958854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
959c5af039cSHong Zhang 
960c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
961c5af039cSHong Zhang   owners = merge->rowmap->range;
962854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
963c5af039cSHong Zhang 
964c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
965c5af039cSHong Zhang   A_loc = ptap->A_loc;
966c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
967c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
968d6ab1dc0SHong Zhang   ai    = a_loc->i;
969d6ab1dc0SHong Zhang   aj    = a_loc->j;
970c5af039cSHong Zhang 
971854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
972c5af039cSHong Zhang 
973c5af039cSHong Zhang   for (i=0; i<am; i++) {
974e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
975d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
976d6ab1dc0SHong Zhang     adj = aj + ai[i];
977d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
978c5af039cSHong Zhang     for (j=0; j<anz; j++) {
979e745005bSHong Zhang       aval[adj[j]] = ada[j];
980c5af039cSHong Zhang     }
981c5af039cSHong Zhang 
982c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
983c5af039cSHong Zhang     /*--------------------------------------------------------------*/
984c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
985c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
986c5af039cSHong Zhang     poJ = po->j + po->i[i];
987c5af039cSHong Zhang     pA  = po->a + po->i[i];
988c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
989c5af039cSHong Zhang       row = poJ[j];
990c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
991c5af039cSHong Zhang       cj  = coj + coi[row];
992c5af039cSHong Zhang       ca  = coa + coi[row];
993c5af039cSHong Zhang       /* perform dense axpy */
994c5af039cSHong Zhang       valtmp = pA[j];
995c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
996e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
997c5af039cSHong Zhang       }
998c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
999c5af039cSHong Zhang     }
1000c5af039cSHong Zhang 
1001c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1002c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1003c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1004c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1005c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1006c5af039cSHong Zhang       row = pdJ[j];
1007c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1008c5af039cSHong Zhang       cj  = bj + bi[row];
1009c5af039cSHong Zhang       ca  = ba + bi[row];
1010c5af039cSHong Zhang       /* perform dense axpy */
1011c5af039cSHong Zhang       valtmp = pA[j];
1012c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1013e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1014c5af039cSHong Zhang       }
1015c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1016c5af039cSHong Zhang     }
1017c5af039cSHong Zhang 
1018d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1019d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1020e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1021c5af039cSHong Zhang   }
1022c5af039cSHong Zhang 
1023c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1024c5af039cSHong Zhang   /*------------------------------------*/
1025c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1026c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1027c5af039cSHong Zhang   len_s  = merge->len_s;
1028c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1029c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1030c5af039cSHong Zhang 
1031dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1032c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1033c5af039cSHong Zhang     if (!len_s[proc]) continue;
1034c5af039cSHong Zhang     i    = merge->owners_co[proc];
1035c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1036c5af039cSHong Zhang     k++;
1037c5af039cSHong Zhang   }
1038c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1039c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1040c5af039cSHong Zhang 
1041c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1042c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1043c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1044c5af039cSHong Zhang 
1045c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1046c5af039cSHong Zhang   /*----------------------------------------------------*/
1047dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1048c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1049c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1050c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1051c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1052c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1053c5af039cSHong Zhang   }
1054c5af039cSHong Zhang 
1055c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1056c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1057c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1058c5af039cSHong Zhang     ba_i = ba + bi[i];
1059c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1060c5af039cSHong Zhang     /* add received vals into ba */
1061c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1062c5af039cSHong Zhang       /* i-th row */
1063c5af039cSHong Zhang       if (i == *nextrow[k]) {
1064c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1065c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1066c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1067c5af039cSHong Zhang         nextcj = 0;
1068c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1069c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1070c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1071c5af039cSHong Zhang           }
1072c5af039cSHong Zhang         }
1073c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1074c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1075c5af039cSHong Zhang       }
1076c5af039cSHong Zhang     }
1077c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1078c5af039cSHong Zhang   }
1079c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1080c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1081c5af039cSHong Zhang 
1082c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1083c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1084c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1085c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1086e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1087f96d37fcSHong Zhang   PetscFunctionReturn(0);
1088f96d37fcSHong Zhang }
1089f96d37fcSHong Zhang 
1090c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1091c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
10922bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1093f96d37fcSHong Zhang {
1094f96d37fcSHong Zhang   PetscErrorCode      ierr;
10954e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1096f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10970298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1098f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1099f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1100f96d37fcSHong Zhang   PetscInt            nnz;
11014e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1102497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1103f96d37fcSHong Zhang   PetscBT             lnkbt;
1104ce94432eSBarry Smith   MPI_Comm            comm;
1105f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1106f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1107f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1108f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1109f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1110f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1111f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1112f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1113d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1114f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1115438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1116f96d37fcSHong Zhang   PetscScalar         *vals;
11174e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1118f96d37fcSHong Zhang 
1119f96d37fcSHong Zhang   PetscFunctionBegin;
1120ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1121c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
11226c4ed002SBarry 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);
1123c5af039cSHong Zhang 
1124f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1125f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1126f96d37fcSHong Zhang 
1127f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1128b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1129f96d37fcSHong Zhang 
1130f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1131f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11322205254eSKarl Rupp 
1133c5af039cSHong Zhang   ptap->A_loc = A_loc;
11342205254eSKarl Rupp 
11351c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1136d6ab1dc0SHong Zhang   ai    = a_loc->i;
1137d6ab1dc0SHong Zhang   aj    = a_loc->j;
1138f96d37fcSHong Zhang 
1139f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1140f96d37fcSHong Zhang   /*----------------------------------------------------*/
11414e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11424e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11434e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11444e938277SHong Zhang 
11454e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11464e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11474e938277SHong Zhang   poti = pot->i; potj = pot->j;
1148f96d37fcSHong Zhang 
1149f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11502205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1151854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1152f96d37fcSHong Zhang   coi[0] = 0;
1153f96d37fcSHong Zhang 
1154f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1155f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1156a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1157f96d37fcSHong Zhang   current_space = free_space;
115819f0e57aSHong Zhang 
1159f96d37fcSHong Zhang   /* create and initialize a linked list */
1160438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1161f96d37fcSHong Zhang 
1162f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1163f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1164f96d37fcSHong Zhang     ptJ = potj + poti[i];
1165f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1166f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1167d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1168d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1169f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1170d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1171f96d37fcSHong Zhang     }
11724e938277SHong Zhang     nnz = lnk[0];
1173f96d37fcSHong Zhang 
1174f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1175f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1176f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1177f96d37fcSHong Zhang       nspacedouble++;
1178f96d37fcSHong Zhang     }
1179f96d37fcSHong Zhang 
1180f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11814e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11822205254eSKarl Rupp 
1183f96d37fcSHong Zhang     current_space->array           += nnz;
1184f96d37fcSHong Zhang     current_space->local_used      += nnz;
1185f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11862205254eSKarl Rupp 
1187f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1188f96d37fcSHong Zhang   }
1189f96d37fcSHong Zhang 
1190854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1191f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11922205254eSKarl Rupp 
1193118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1194f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1195f96d37fcSHong Zhang 
1196f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1197f96d37fcSHong Zhang   /*----------------------------------------------*/
1198f96d37fcSHong Zhang   /* determine row ownership */
1199b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1200f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12012205254eSKarl Rupp 
1202f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1203f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12042205254eSKarl Rupp 
1205f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1206f96d37fcSHong Zhang   owners = merge->rowmap->range;
1207f96d37fcSHong Zhang 
1208f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12091795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1210785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12112205254eSKarl Rupp 
1212f96d37fcSHong Zhang   len_s        = merge->len_s;
1213f96d37fcSHong Zhang   merge->nsend = 0;
1214f96d37fcSHong Zhang 
1215854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1216f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1217f96d37fcSHong Zhang 
1218f96d37fcSHong Zhang   proc = 0;
1219f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1220f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1221f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1222f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1223f96d37fcSHong Zhang   }
1224f96d37fcSHong Zhang 
1225f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1226f96d37fcSHong Zhang   owners_co[0] = 0;
1227f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1228f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1229f96d37fcSHong Zhang     if (len_si[proc]) {
1230f96d37fcSHong Zhang       merge->nsend++;
1231f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1232f96d37fcSHong Zhang       len         += len_si[proc];
1233f96d37fcSHong Zhang     }
1234f96d37fcSHong Zhang   }
1235f96d37fcSHong Zhang 
1236f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12370298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1238f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1239f96d37fcSHong Zhang 
1240f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1241f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1242f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1243854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1244f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1245f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1246f96d37fcSHong Zhang     i    = owners_co[proc];
1247f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1248f96d37fcSHong Zhang     k++;
1249f96d37fcSHong Zhang   }
1250f96d37fcSHong Zhang 
1251f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1252785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1253f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1254c280ed6aSJed Brown     PetscMPIInt icompleted;
1255c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1256f96d37fcSHong Zhang   }
1257f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1258f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1259f96d37fcSHong Zhang 
1260f96d37fcSHong Zhang   /* send and recv coi */
1261f96d37fcSHong Zhang   /*-------------------*/
1262f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1263f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1264854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1265f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1266f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1267f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1268f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1269f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1270f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1271f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1272f96d37fcSHong Zhang     */
1273f96d37fcSHong Zhang     /*-------------------------------------------*/
1274f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1275f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1276f96d37fcSHong Zhang     buf_si[0]   = nrows;
1277f96d37fcSHong Zhang     buf_si_i[0] = 0;
1278f96d37fcSHong Zhang     nrows       = 0;
1279f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1280f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1281f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1282f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1283f96d37fcSHong Zhang       nrows++;
1284f96d37fcSHong Zhang     }
1285f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1286f96d37fcSHong Zhang     k++;
1287f96d37fcSHong Zhang     buf_si += len_si[proc];
1288f96d37fcSHong Zhang   }
1289f96d37fcSHong Zhang   i = merge->nrecv;
1290f96d37fcSHong Zhang   while (i--) {
1291c280ed6aSJed Brown     PetscMPIInt icompleted;
1292c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1293f96d37fcSHong Zhang   }
1294f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1295f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1296f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1297f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1298f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1299f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1300f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1301f96d37fcSHong Zhang 
1302f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1303f96d37fcSHong Zhang   /*------------------------------------------*/
1304f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1305854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1306f96d37fcSHong Zhang   bi[0] = 0;
1307f96d37fcSHong Zhang 
1308c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1309f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1310a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1311f96d37fcSHong Zhang   current_space = free_space;
1312f96d37fcSHong Zhang 
1313dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1314f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1315f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1316f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1317f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1318f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1319f96d37fcSHong Zhang   }
1320f96d37fcSHong Zhang 
13211c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1322f96d37fcSHong Zhang   rmax = 0;
1323f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1324f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1325f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1326f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1327f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1328f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1329d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1330d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1331f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1332d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1333f96d37fcSHong Zhang     }
13344e938277SHong Zhang 
1335f96d37fcSHong Zhang     /* add received col data into lnk */
1336f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1337f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1338f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1339f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13404e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1341f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1342f96d37fcSHong Zhang       }
1343f96d37fcSHong Zhang     }
13444e938277SHong Zhang     nnz = lnk[0];
1345f96d37fcSHong Zhang 
1346f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1347f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1348f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1349f96d37fcSHong Zhang       nspacedouble++;
1350f96d37fcSHong Zhang     }
1351f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13524e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1353f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13542205254eSKarl Rupp 
1355f96d37fcSHong Zhang     current_space->array           += nnz;
1356f96d37fcSHong Zhang     current_space->local_used      += nnz;
1357f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13582205254eSKarl Rupp 
1359f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1360f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1361f96d37fcSHong Zhang   }
1362f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1363f96d37fcSHong Zhang 
1364854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1365f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13662205254eSKarl Rupp 
1367118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1368f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1369d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13704e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13714e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1372f96d37fcSHong Zhang 
13731c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13741c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13751795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1376f96d37fcSHong Zhang 
1377f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1378f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
137933d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1380f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1381f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1382f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1383f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1384f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1385f96d37fcSHong Zhang     row  = i + rstart;
1386f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1387f96d37fcSHong Zhang     Jptr = bj + bi[i];
1388f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1389f96d37fcSHong Zhang   }
1390f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1391f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1392f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1393f96d37fcSHong Zhang 
1394f96d37fcSHong Zhang   merge->bi        = bi;
1395f96d37fcSHong Zhang   merge->bj        = bj;
1396f96d37fcSHong Zhang   merge->coi       = coi;
1397f96d37fcSHong Zhang   merge->coj       = coj;
1398f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1399f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1400f96d37fcSHong Zhang   merge->owners_co = owners_co;
1401f96d37fcSHong Zhang 
1402f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1403f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1404f96d37fcSHong Zhang   c->ptap     = ptap;
14050298fd71SBarry Smith   ptap->api   = NULL;
14060298fd71SBarry Smith   ptap->apj   = NULL;
1407f96d37fcSHong Zhang   ptap->merge = merge;
1408375ed354SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1409375ed354SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1410375ed354SHong Zhang 
1411375ed354SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1412375ed354SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1413375ed354SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1414d6ab1dc0SHong Zhang 
1415d6ab1dc0SHong Zhang   *C = Cmpi;
1416d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1417d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
141857622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
141957622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1420d6ab1dc0SHong Zhang   } else {
1421d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1422d6ab1dc0SHong Zhang   }
1423d6ab1dc0SHong Zhang #endif
1424d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1425d6ab1dc0SHong Zhang }
1426d6ab1dc0SHong Zhang 
14276da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1428d6ab1dc0SHong Zhang {
1429d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1430d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1431d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1432d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1433d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1434e745005bSHong Zhang   PetscInt            *adj;
1435e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1436e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1437d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1438ce94432eSBarry Smith   MPI_Comm            comm;
1439d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1440d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1441d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1442d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1443d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1444d6ab1dc0SHong Zhang   MPI_Status          *status;
1445d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1446d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1447a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1448d6ab1dc0SHong Zhang   Mat                 A_loc;
1449d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1450d6ab1dc0SHong Zhang 
1451d6ab1dc0SHong Zhang   PetscFunctionBegin;
1452ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1453d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1454d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1455d6ab1dc0SHong Zhang 
1456d6ab1dc0SHong Zhang   ptap  = c->ptap;
1457d6ab1dc0SHong Zhang   merge = ptap->merge;
1458d6ab1dc0SHong Zhang 
1459e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1460e745005bSHong Zhang   /*------------------------------------------*/
1461d6ab1dc0SHong Zhang   /* get data from symbolic products */
1462d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1463854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1464d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1465d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14661795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1467d6ab1dc0SHong Zhang 
1468d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1469d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1470d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1471d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1472d6ab1dc0SHong Zhang   ai    = a_loc->i;
1473d6ab1dc0SHong Zhang   aj    = a_loc->j;
1474d6ab1dc0SHong Zhang 
1475d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1476d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1477d6ab1dc0SHong Zhang     adj = aj + ai[i];
1478d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1479d6ab1dc0SHong Zhang 
1480d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1481e745005bSHong Zhang     /*-------------------------------------------------------------*/
1482d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1483d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1484d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1485d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1486d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1487d6ab1dc0SHong Zhang       row = poJ[j];
1488d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1489d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1490e745005bSHong Zhang       /* perform sparse axpy */
1491e745005bSHong Zhang       nexta  = 0;
1492d6ab1dc0SHong Zhang       valtmp = pA[j];
1493e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1494e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1495e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1496e745005bSHong Zhang           nexta++;
1497d6ab1dc0SHong Zhang         }
1498e745005bSHong Zhang       }
1499e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1500d6ab1dc0SHong Zhang     }
1501d6ab1dc0SHong Zhang 
1502d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1503d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1504d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1505d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1506d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1507d6ab1dc0SHong Zhang       row = pdJ[j];
1508d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1509d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1510e745005bSHong Zhang       /* perform sparse axpy */
1511e745005bSHong Zhang       nexta  = 0;
1512d6ab1dc0SHong Zhang       valtmp = pA[j];
1513e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1514e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1515e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1516e745005bSHong Zhang           nexta++;
1517d6ab1dc0SHong Zhang         }
1518e745005bSHong Zhang       }
1519e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1520d6ab1dc0SHong Zhang     }
1521d6ab1dc0SHong Zhang   }
1522d6ab1dc0SHong Zhang 
1523d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1524d6ab1dc0SHong Zhang   /*------------------------------------*/
1525d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1526d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1527d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1528d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1530d6ab1dc0SHong Zhang 
1531dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1532d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1533d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1534d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1535d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1536d6ab1dc0SHong Zhang     k++;
1537d6ab1dc0SHong Zhang   }
1538d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1539d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1540d6ab1dc0SHong Zhang 
1541d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1542d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1543d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1544d6ab1dc0SHong Zhang 
1545d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1546d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1547dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1548d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1549e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1550d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1551d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1552d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1553d6ab1dc0SHong Zhang   }
1554d6ab1dc0SHong Zhang 
1555d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1556d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1557d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1558d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1559d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1560d6ab1dc0SHong Zhang     /* add received vals into ba */
1561d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1562d6ab1dc0SHong Zhang       /* i-th row */
1563d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1564d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1565d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1566d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1567d6ab1dc0SHong Zhang         nextcj = 0;
1568d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1569d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1570d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1571d6ab1dc0SHong Zhang           }
1572d6ab1dc0SHong Zhang         }
1573d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1574d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1575d6ab1dc0SHong Zhang       }
1576d6ab1dc0SHong Zhang     }
1577d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1578d6ab1dc0SHong Zhang   }
1579d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1580d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1581d6ab1dc0SHong Zhang 
1582d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1583d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1584d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1585d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1586d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1587d6ab1dc0SHong Zhang }
1588d6ab1dc0SHong Zhang 
1589c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
15902bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15912bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
15926da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1593d6ab1dc0SHong Zhang {
1594d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1595d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1596d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15970298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
159833ba5e2fSHong Zhang   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1599d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1600d6ab1dc0SHong Zhang   PetscInt            nnz;
1601d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1602d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1603ce94432eSBarry Smith   MPI_Comm            comm;
1604d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1605d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1606d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1607d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1608d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1609d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1610d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1611d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1612d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1613d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
16144b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1615d6ab1dc0SHong Zhang   PetscScalar         *vals;
1616d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
16170acc65b4SHong Zhang   PetscTable          ta;
1618d6ab1dc0SHong Zhang 
1619d6ab1dc0SHong Zhang   PetscFunctionBegin;
1620ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1621d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
16226c4ed002SBarry 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);
1623d6ab1dc0SHong Zhang 
1624d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1625d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1626d6ab1dc0SHong Zhang 
1627d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1628b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1629d6ab1dc0SHong Zhang 
1630d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1631d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16322205254eSKarl Rupp 
1633d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1634d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1635d6ab1dc0SHong Zhang   ai          = a_loc->i;
1636d6ab1dc0SHong Zhang   aj          = a_loc->j;
1637d6ab1dc0SHong Zhang 
1638d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1639d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1640d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1641d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1642d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1643d6ab1dc0SHong Zhang 
1644d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1645d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1646d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1647d6ab1dc0SHong Zhang 
1648d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1649d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1650d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1651854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1652d6ab1dc0SHong Zhang   coi[0] = 0;
1653d6ab1dc0SHong Zhang 
1654d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1655f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1656a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1657d6ab1dc0SHong Zhang   current_space = free_space;
165819f0e57aSHong Zhang 
1659d6ab1dc0SHong Zhang   /* create and initialize a linked list */
166033ba5e2fSHong Zhang   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
16614b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
16620acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
166333ba5e2fSHong Zhang 
16640acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1665d6ab1dc0SHong Zhang 
1666d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1667d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1668d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1669d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1670d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1671d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1672d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1673d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1674d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1675d6ab1dc0SHong Zhang     }
1676d6ab1dc0SHong Zhang     nnz = lnk[0];
1677d6ab1dc0SHong Zhang 
1678d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1679d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1680f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1681d6ab1dc0SHong Zhang       nspacedouble++;
1682d6ab1dc0SHong Zhang     }
1683d6ab1dc0SHong Zhang 
1684d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1685d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16862205254eSKarl Rupp 
1687d6ab1dc0SHong Zhang     current_space->array           += nnz;
1688d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1689d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16902205254eSKarl Rupp 
1691d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1692d6ab1dc0SHong Zhang   }
1693d6ab1dc0SHong Zhang 
1694854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1695d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16960acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
16972205254eSKarl Rupp 
1698118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1699d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1700d6ab1dc0SHong Zhang 
1701d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1702d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1703d6ab1dc0SHong Zhang   /* determine row ownership */
1704b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1705d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17062205254eSKarl Rupp 
1707d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1708d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17092205254eSKarl Rupp 
1710d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1711d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1712d6ab1dc0SHong Zhang 
1713d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17141795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1715785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17162205254eSKarl Rupp 
1717d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1718d6ab1dc0SHong Zhang   merge->nsend = 0;
1719d6ab1dc0SHong Zhang 
1720854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1721d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1722d6ab1dc0SHong Zhang 
1723d6ab1dc0SHong Zhang   proc = 0;
1724d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1725d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1726d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1727d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1728d6ab1dc0SHong Zhang   }
1729d6ab1dc0SHong Zhang 
1730d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1731d6ab1dc0SHong Zhang   owners_co[0] = 0;
1732d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1733d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1734d6ab1dc0SHong Zhang     if (len_si[proc]) {
1735d6ab1dc0SHong Zhang       merge->nsend++;
1736d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1737d6ab1dc0SHong Zhang       len         += len_si[proc];
1738d6ab1dc0SHong Zhang     }
1739d6ab1dc0SHong Zhang   }
1740d6ab1dc0SHong Zhang 
1741d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17420298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1743d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang 
1745d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1746d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1747d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1748854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1749d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1750d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1751d6ab1dc0SHong Zhang     i    = owners_co[proc];
1752d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1753d6ab1dc0SHong Zhang     k++;
1754d6ab1dc0SHong Zhang   }
1755d6ab1dc0SHong Zhang 
1756d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1757785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1758d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1759c280ed6aSJed Brown     PetscMPIInt icompleted;
1760c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   }
1762d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1764d6ab1dc0SHong Zhang 
17650acc65b4SHong Zhang   /* add received column indices into table to update Armax */
176633ba5e2fSHong Zhang   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
17670acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
17680acc65b4SHong Zhang     Jptr = buf_rj[k];
17690acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
17700acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
17710acc65b4SHong Zhang     }
17720acc65b4SHong Zhang   }
17730acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
177433ba5e2fSHong 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); */
17750acc65b4SHong Zhang 
1776d6ab1dc0SHong Zhang   /* send and recv coi */
1777d6ab1dc0SHong Zhang   /*-------------------*/
1778d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1779d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1780854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1781d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1782d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1783d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1784d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1785d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1786d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1787d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1788d6ab1dc0SHong Zhang     */
1789d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1790d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1791d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1792d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1793d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1794d6ab1dc0SHong Zhang     nrows       = 0;
1795d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1796d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1797d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1798d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1799d6ab1dc0SHong Zhang       nrows++;
1800d6ab1dc0SHong Zhang     }
1801d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1802d6ab1dc0SHong Zhang     k++;
1803d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1804d6ab1dc0SHong Zhang   }
1805d6ab1dc0SHong Zhang   i = merge->nrecv;
1806d6ab1dc0SHong Zhang   while (i--) {
1807c280ed6aSJed Brown     PetscMPIInt icompleted;
1808c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1809d6ab1dc0SHong Zhang   }
1810d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1811d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1812d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1813d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1815d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1816d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1817d6ab1dc0SHong Zhang 
1818d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1819d6ab1dc0SHong Zhang   /*------------------------------------------*/
1820d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1821854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1822d6ab1dc0SHong Zhang   bi[0] = 0;
1823d6ab1dc0SHong Zhang 
1824d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1825f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1826a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1827d6ab1dc0SHong Zhang   current_space = free_space;
1828d6ab1dc0SHong Zhang 
1829dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1830d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1831d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1832d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1833d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18342205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1835d6ab1dc0SHong Zhang   }
1836d6ab1dc0SHong Zhang 
18370acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1838d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1839d6ab1dc0SHong Zhang   rmax = 0;
1840d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1841d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1842d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1843d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1844d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1845d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1846d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1847d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1848d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1849d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1850d6ab1dc0SHong Zhang     }
1851d6ab1dc0SHong Zhang 
1852d6ab1dc0SHong Zhang     /* add received col data into lnk */
1853d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1854d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1855d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1856d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1857d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1858d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1859d6ab1dc0SHong Zhang       }
1860d6ab1dc0SHong Zhang     }
1861d6ab1dc0SHong Zhang     nnz = lnk[0];
1862d6ab1dc0SHong Zhang 
1863d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1864d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1865f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1866d6ab1dc0SHong Zhang       nspacedouble++;
1867d6ab1dc0SHong Zhang     }
1868d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1869d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1870d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18712205254eSKarl Rupp 
1872d6ab1dc0SHong Zhang     current_space->array           += nnz;
1873d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1874d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18752205254eSKarl Rupp 
1876d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1877d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1878d6ab1dc0SHong Zhang   }
1879f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang 
1881854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1882d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1883118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1884d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1885d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
18860acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
18870acc65b4SHong Zhang 
1888d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1889d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1890d6ab1dc0SHong Zhang 
1891d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1892d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18931795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1894d6ab1dc0SHong Zhang 
1895d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1896d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
189733d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1898d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1899d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1900d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1901d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1902d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1903d6ab1dc0SHong Zhang     row  = i + rstart;
1904d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1905d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1906d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1907d6ab1dc0SHong Zhang   }
1908d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1909d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1910d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1911d6ab1dc0SHong Zhang 
1912d6ab1dc0SHong Zhang   merge->bi        = bi;
1913d6ab1dc0SHong Zhang   merge->bj        = bj;
1914d6ab1dc0SHong Zhang   merge->coi       = coi;
1915d6ab1dc0SHong Zhang   merge->coj       = coj;
1916d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1917d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1918d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1919d6ab1dc0SHong Zhang 
1920d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1921d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19222205254eSKarl Rupp 
1923d6ab1dc0SHong Zhang   c->ptap     = ptap;
19240298fd71SBarry Smith   ptap->api   = NULL;
19250298fd71SBarry Smith   ptap->apj   = NULL;
1926d6ab1dc0SHong Zhang   ptap->merge = merge;
19270298fd71SBarry Smith   ptap->apa   = NULL;
1928a560ef98SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1929a560ef98SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1930a560ef98SHong Zhang 
1931a560ef98SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1932a560ef98SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1933a560ef98SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1934f96d37fcSHong Zhang 
1935f96d37fcSHong Zhang   *C = Cmpi;
1936f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1937f96d37fcSHong Zhang   if (bi[pn] != 0) {
193857622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
193957622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1940f96d37fcSHong Zhang   } else {
1941f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1942f96d37fcSHong Zhang   }
1943f96d37fcSHong Zhang #endif
1944f96d37fcSHong Zhang   PetscFunctionReturn(0);
1945f96d37fcSHong Zhang }
1946