xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 6c4ed00291fe44f94936dd2f04c02ab3c442e77c)
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 
132499ec78SHong Zhang #undef __FUNCT__
142499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
152499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
162499ec78SHong Zhang {
172499ec78SHong Zhang   PetscErrorCode ierr;
180fc8cf34SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
190d3441aeSHong Zhang   PetscInt       alg=1; /* set default algorithm */
20e55be12dSHong Zhang   MPI_Comm       comm;
212499ec78SHong Zhang 
222499ec78SHong Zhang   PetscFunctionBegin;
232499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
24e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
25*6c4ed002SBarry 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);
26e55be12dSHong Zhang 
270fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
2878d30b94SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr);
290fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
300fc8cf34SHong Zhang 
313ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
320fc8cf34SHong Zhang     switch (alg) {
330fc8cf34SHong Zhang     case 1:
340fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
350fc8cf34SHong Zhang       break;
360fc8cf34SHong Zhang     default:
37b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
380fc8cf34SHong Zhang       break;
390fc8cf34SHong Zhang     }
403ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
412499ec78SHong Zhang   }
423ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
43598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
443ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
452499ec78SHong Zhang   PetscFunctionReturn(0);
462499ec78SHong Zhang }
472499ec78SHong Zhang 
48f33d1a9aSHong Zhang #undef __FUNCT__
49a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
50a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
51598bc09dSHong Zhang {
52598bc09dSHong Zhang   PetscErrorCode ierr;
53598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
54598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
55598bc09dSHong Zhang 
56598bc09dSHong Zhang   PetscFunctionBegin;
57b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
58598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
60a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
61ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
62a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
63a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
64598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
65598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
66598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
67598bc09dSHong Zhang   PetscFunctionReturn(0);
68598bc09dSHong Zhang }
69598bc09dSHong Zhang 
702499ec78SHong Zhang #undef __FUNCT__
71a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
72a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
734ae0847bSHong Zhang {
744ae0847bSHong Zhang   PetscErrorCode ierr;
754ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
764ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
774ae0847bSHong Zhang 
784ae0847bSHong Zhang   PetscFunctionBegin;
794ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
8026fbe8dcSKarl Rupp 
814ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
824ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
834ae0847bSHong Zhang   PetscFunctionReturn(0);
844ae0847bSHong Zhang }
854ae0847bSHong Zhang 
864ae0847bSHong Zhang #undef __FUNCT__
87f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
88f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
89598bc09dSHong Zhang {
90598bc09dSHong Zhang   PetscErrorCode ierr;
914ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
92598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
93598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
94c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
95598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
969ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
97c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
98598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
99c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
10029825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
101598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1029467f45fSHong Zhang   MPI_Comm       comm;
1039467f45fSHong Zhang   PetscMPIInt    size;
104598bc09dSHong Zhang 
105598bc09dSHong Zhang   PetscFunctionBegin;
1069467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1079467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1089467f45fSHong Zhang 
109a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
110598bc09dSHong Zhang   /*-----------------------------------------------------*/
111a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
112b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
113a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
114598bc09dSHong Zhang 
115598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
116598bc09dSHong Zhang   /*----------------------------------------------------------*/
117598bc09dSHong Zhang   /* get data from symbolic products */
118a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
119c12557a1SHong Zhang   p_oth = NULL;
1209467f45fSHong Zhang   if (size >1) {
1219467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1229467f45fSHong Zhang   }
123598bc09dSHong Zhang 
124598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
125598bc09dSHong Zhang   apa = ptap->apa;
126598bc09dSHong Zhang 
12729825010SHong Zhang   api = ptap->api;
12829825010SHong Zhang   apj = ptap->apj;
129598bc09dSHong Zhang   for (i=0; i<cm; i++) {
130c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
131e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
132598bc09dSHong Zhang 
133598bc09dSHong Zhang     /* set values in C */
134598bc09dSHong Zhang     apJ  = apj + api[i];
135598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
136598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
137598bc09dSHong Zhang 
138598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
139598bc09dSHong Zhang     ca = coa + co->i[i];
140598bc09dSHong Zhang     k  = 0;
141598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
142598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
143598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
144598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
145598bc09dSHong Zhang       k++;
146598bc09dSHong Zhang     }
147598bc09dSHong Zhang 
148598bc09dSHong Zhang     /* diagonal part of C */
149598bc09dSHong Zhang     ca = cda + cd->i[i];
150598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
151598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
152598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
153598bc09dSHong Zhang       k++;
154598bc09dSHong Zhang     }
155598bc09dSHong Zhang 
156598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
157598bc09dSHong Zhang     ca = coa + co->i[i];
158598bc09dSHong Zhang     for (; k0<conz; k0++) {
159598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
160598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
161598bc09dSHong Zhang       k++;
162598bc09dSHong Zhang     }
163598bc09dSHong Zhang   }
164598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166598bc09dSHong Zhang   PetscFunctionReturn(0);
167598bc09dSHong Zhang }
168598bc09dSHong Zhang 
169598bc09dSHong Zhang #undef __FUNCT__
1700fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1710fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
172598bc09dSHong Zhang {
173598bc09dSHong Zhang   PetscErrorCode     ierr;
174ce94432eSBarry Smith   MPI_Comm           comm;
1759467f45fSHong Zhang   PetscMPIInt        size;
176bfcd1a73SHong Zhang   Mat                Cmpi;
177598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1780298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1794ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
180bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1814ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1824ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
183d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
184e55be12dSHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n,Crmax;
185598bc09dSHong Zhang   PetscBT            lnkbt;
186598bc09dSHong Zhang   PetscScalar        *apa;
187bfcd1a73SHong Zhang   PetscReal          afill;
188e55be12dSHong Zhang   PetscTable         ta;
189598bc09dSHong Zhang 
190598bc09dSHong Zhang   PetscFunctionBegin;
191ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1929467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1939467f45fSHong Zhang 
194598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
195b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
196598bc09dSHong Zhang 
197598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
198b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
19919f0e57aSHong Zhang 
200598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
201a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
202598bc09dSHong Zhang 
203a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
204598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2059467f45fSHong Zhang   if (size > 1) {
2069467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
207598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
20820e3dc75SHong Zhang   } else {
20920e3dc75SHong Zhang     p_oth = NULL;
21020e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2119467f45fSHong Zhang   }
212598bc09dSHong Zhang 
213598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
214598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
215854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
216a1a4e84aSHong Zhang   ptap->api = api;
217598bc09dSHong Zhang   api[0]    = 0;
218598bc09dSHong Zhang 
219598bc09dSHong Zhang   /* create and initialize a linked list */
220e717af80SHong Zhang   Crmax = 6*(p_loc->rmax + (PetscInt)(1.e-2*pN));
221eca6b86aSHong Zhang   if (Crmax > pN) Crmax = pN;
222eca6b86aSHong Zhang   ierr = PetscTableCreate(Crmax,pN,&ta);CHKERRQ(ierr);
2234b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
2244b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
225e55be12dSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
226e55be12dSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
227e55be12dSHong Zhang 
228e55be12dSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
229598bc09dSHong Zhang 
230bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
231bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
232598bc09dSHong Zhang   current_space = free_space;
233598bc09dSHong Zhang 
234598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
235598bc09dSHong Zhang   for (i=0; i<am; i++) {
236598bc09dSHong Zhang     /* diagonal portion of A */
237598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
238598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
239598bc09dSHong Zhang       row  = *adj++;
240598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
241598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
242598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2431fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
244598bc09dSHong Zhang     }
245598bc09dSHong Zhang     /* off-diagonal portion of A */
246598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
247598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
248598bc09dSHong Zhang       row  = *aoj++;
249598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
250598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2511fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
252598bc09dSHong Zhang     }
253598bc09dSHong Zhang 
254d14fa243SHong Zhang     apnz     = lnk[0];
255598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
256598bc09dSHong Zhang 
257598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
258598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
259598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
260598bc09dSHong Zhang       nspacedouble++;
261598bc09dSHong Zhang     }
262598bc09dSHong Zhang 
263598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2641fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
265598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2662205254eSKarl Rupp 
267598bc09dSHong Zhang     current_space->array           += apnz;
268598bc09dSHong Zhang     current_space->local_used      += apnz;
269598bc09dSHong Zhang     current_space->local_remaining -= apnz;
270598bc09dSHong Zhang   }
271598bc09dSHong Zhang 
272598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
273598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
274854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
275a1a4e84aSHong Zhang   apj  = ptap->apj;
276a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
277598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
278598bc09dSHong Zhang 
279f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2801795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2812205254eSKarl Rupp 
282f84c536bSHong Zhang   ptap->apa = apa;
283f84c536bSHong Zhang 
284bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
285598bc09dSHong Zhang   /*----------------------------------------------------*/
286bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
287bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
28833d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
289a2f3521dSMark F. Adams 
290bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
291bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
292598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
293598bc09dSHong Zhang   for (i=0; i<am; i++) {
294598bc09dSHong Zhang     row  = i + rstart;
295598bc09dSHong Zhang     apnz = api[i+1] - api[i];
296bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
297598bc09dSHong Zhang     apj += apnz;
298598bc09dSHong Zhang   }
299bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
300bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
301598bc09dSHong Zhang 
302bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
303bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
304f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
305bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
306bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
307598bc09dSHong Zhang 
308bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
309bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
310598bc09dSHong Zhang   c->ptap = ptap;
311598bc09dSHong Zhang 
312bfcd1a73SHong Zhang   *C = Cmpi;
313bfcd1a73SHong Zhang 
314bfcd1a73SHong Zhang   /* set MatInfo */
315118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
316bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
317bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
318bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
319bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
320bfcd1a73SHong Zhang 
321bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
322bfcd1a73SHong Zhang   if (api[am]) {
32357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
32457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
325bfcd1a73SHong Zhang   } else {
326bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
327bfcd1a73SHong Zhang   }
328bfcd1a73SHong Zhang #endif
329598bc09dSHong Zhang   PetscFunctionReturn(0);
330598bc09dSHong Zhang }
331598bc09dSHong Zhang 
3329123193aSHong Zhang #undef __FUNCT__
3339123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3349123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3359123193aSHong Zhang {
3369123193aSHong Zhang   PetscErrorCode ierr;
3379123193aSHong Zhang 
3389123193aSHong Zhang   PetscFunctionBegin;
3399123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3403ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3419123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3423ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3439123193aSHong Zhang   }
3443ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3459123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3463ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3479123193aSHong Zhang   PetscFunctionReturn(0);
3489123193aSHong Zhang }
3499123193aSHong Zhang 
3504324174eSBarry Smith typedef struct {
3514324174eSBarry Smith   Mat         workB;
3524324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3534324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3544324174eSBarry Smith } MPIAIJ_MPIDense;
3554324174eSBarry Smith 
3567af9e4fdSHong Zhang #undef __FUNCT__
35796b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
35896b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3594324174eSBarry Smith {
3604324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3614324174eSBarry Smith   PetscErrorCode  ierr;
3624324174eSBarry Smith 
3634324174eSBarry Smith   PetscFunctionBegin;
3646bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3654324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3664324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3674324174eSBarry Smith   PetscFunctionReturn(0);
3684324174eSBarry Smith }
3694324174eSBarry Smith 
3709123193aSHong Zhang #undef __FUNCT__
371fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense"
372fd4e9aacSBarry Smith /*
373fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
374fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
375fd4e9aacSBarry Smith 
376fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
377fd4e9aacSBarry Smith */
378fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
379fd4e9aacSBarry Smith {
380fd4e9aacSBarry Smith   PetscErrorCode         ierr;
381fd4e9aacSBarry Smith   PetscBool              flg;
382fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
383fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
384fd4e9aacSBarry Smith   PetscContainer         container;
385fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
386fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
387fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
388fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
389fd4e9aacSBarry Smith 
390fd4e9aacSBarry Smith   PetscFunctionBegin;
391fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
392fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
393fd4e9aacSBarry Smith 
394fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
395fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
396fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
397fd4e9aacSBarry Smith 
398fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
399fd4e9aacSBarry Smith 
400fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
401fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
402fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
403fd4e9aacSBarry Smith   /* Create work arrays needed */
404fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
405fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
406fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
407fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
408fd4e9aacSBarry Smith 
409fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
410fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
411fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
412fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
413fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
414fd4e9aacSBarry Smith 
415fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
416fd4e9aacSBarry Smith   PetscFunctionReturn(0);
417fd4e9aacSBarry Smith }
418fd4e9aacSBarry Smith 
419fd4e9aacSBarry Smith #undef __FUNCT__
4209123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4219123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4229123193aSHong Zhang {
4239123193aSHong Zhang   PetscErrorCode         ierr;
424f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
425d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
426bf0cc555SLisandro Dalcin   PetscContainer         container;
4274324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4284324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4294324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4304324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
431d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4329123193aSHong Zhang 
4339123193aSHong Zhang   PetscFunctionBegin;
434ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
435d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
43633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
437cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4380298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
439cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
440cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4412205254eSKarl Rupp 
442d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
443f32d5d43SBarry Smith 
444b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
445f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4460298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4474324174eSBarry Smith   /* Create work arrays needed */
448dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
449dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
450dcca6d9dSJed Brown                       from->n,&contents->rwaits,
451dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4524324174eSBarry Smith 
453ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
454bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
45596b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
456bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
457bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4589123193aSHong Zhang   PetscFunctionReturn(0);
4599123193aSHong Zhang }
4609123193aSHong Zhang 
4617af9e4fdSHong Zhang #undef __FUNCT__
4627af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
463f32d5d43SBarry Smith /*
4642636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4652636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
466f32d5d43SBarry Smith */
4674324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
468f32d5d43SBarry Smith {
469f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
470f32d5d43SBarry Smith   PetscErrorCode         ierr;
471f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
472f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
473f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
474f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
475f32d5d43SBarry Smith   PetscInt               i,j,k;
476aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
477aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
478f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
479ce94432eSBarry Smith   MPI_Comm               comm;
480d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
481f32d5d43SBarry Smith   MPI_Status             status;
4824324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
483bf0cc555SLisandro Dalcin   PetscContainer         container;
4844324174eSBarry Smith   Mat                    workB;
485f32d5d43SBarry Smith 
486f32d5d43SBarry Smith   PetscFunctionBegin;
487ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
488bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
48929825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
490bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4914324174eSBarry Smith 
4924324174eSBarry Smith   workB = *outworkB = contents->workB;
493cf1a0672SHong 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);
494f32d5d43SBarry Smith   sindices = to->indices;
495f32d5d43SBarry Smith   sstarts  = to->starts;
496f32d5d43SBarry Smith   sprocs   = to->procs;
4974324174eSBarry Smith   swaits   = contents->swaits;
4984324174eSBarry Smith   svalues  = contents->svalues;
499f32d5d43SBarry Smith 
500f32d5d43SBarry Smith   rindices = from->indices;
501f32d5d43SBarry Smith   rstarts  = from->starts;
502f32d5d43SBarry Smith   rprocs   = from->procs;
5034324174eSBarry Smith   rwaits   = contents->rwaits;
5044324174eSBarry Smith   rvalues  = contents->rvalues;
505f32d5d43SBarry Smith 
5068c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5078c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
508f32d5d43SBarry Smith 
509f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
510f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
511f32d5d43SBarry Smith   }
512f32d5d43SBarry Smith 
513f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
514f32d5d43SBarry Smith     /* pack a message at a time */
515f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
516f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
517f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5182636ff26SBarry Smith       }
5192636ff26SBarry Smith     }
520f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
521f32d5d43SBarry Smith   }
5222636ff26SBarry Smith 
523f32d5d43SBarry Smith   nrecvs = from->n;
524f32d5d43SBarry Smith   while (nrecvs) {
525f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
526f32d5d43SBarry Smith     nrecvs--;
527f32d5d43SBarry Smith     /* unpack a message at a time */
528f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
529f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
530f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5312636ff26SBarry Smith       }
5322636ff26SBarry Smith     }
533f32d5d43SBarry Smith   }
534cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
535f32d5d43SBarry Smith 
5368c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5378c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
538f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
539f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
540f32d5d43SBarry Smith   PetscFunctionReturn(0);
541f32d5d43SBarry Smith }
542f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
543f32d5d43SBarry Smith 
5449123193aSHong Zhang #undef __FUNCT__
5459123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
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 
5681634b032SHong Zhang #undef __FUNCT__
569f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
570f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5711634b032SHong Zhang {
5721634b032SHong Zhang   PetscErrorCode ierr;
573cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
574cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
575cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
576cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
577cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
578cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
579cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
58029825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
58129825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
582cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
58329825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
584cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
58529825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
58629825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
587a7c7454dSHong Zhang   MPI_Comm       comm;
588a7c7454dSHong Zhang   PetscMPIInt    size;
5891634b032SHong Zhang 
5901634b032SHong Zhang   PetscFunctionBegin;
591a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
592a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
593a7c7454dSHong Zhang 
594cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
595cf1a0672SHong Zhang   /*-----------------------------------------------------*/
596cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
597b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
598cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
599cf1a0672SHong Zhang 
600cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
601cf1a0672SHong Zhang   /*----------------------------------------------------------*/
602cf1a0672SHong Zhang   /* get data from symbolic products */
603cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
604cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
605a7c7454dSHong Zhang   if (size >1) {
606a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
607cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
60820e3dc75SHong Zhang   } else {
60920e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
610a7c7454dSHong Zhang   }
611cf1a0672SHong Zhang 
61229825010SHong Zhang   api = ptap->api;
61329825010SHong Zhang   apj = ptap->apj;
614cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
61529825010SHong Zhang     apJ = apj + api[i];
61629825010SHong Zhang 
617cf1a0672SHong Zhang     /* diagonal portion of A */
618cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
619cf1a0672SHong Zhang     adj = ad->j + adi[i];
620cf1a0672SHong Zhang     ada = ad->a + adi[i];
621cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
622cf1a0672SHong Zhang       row = adj[j];
623cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
624cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
625cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
62629825010SHong Zhang       /* perform sparse axpy */
627cf1a0672SHong Zhang       valtmp = ada[j];
62829825010SHong Zhang       nextp  = 0;
62929825010SHong Zhang       for (k=0; nextp<pnz; k++) {
63029825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
63129825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63229825010SHong Zhang         }
6331634b032SHong Zhang       }
634cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
635cf1a0672SHong Zhang     }
6361634b032SHong Zhang 
637cf1a0672SHong Zhang     /* off-diagonal portion of A */
638cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
639cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
640cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
641cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
642cf1a0672SHong Zhang       row = aoj[j];
643cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
644cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
645cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
64629825010SHong Zhang       /* perform sparse axpy */
647cf1a0672SHong Zhang       valtmp = aoa[j];
64829825010SHong Zhang       nextp  = 0;
64929825010SHong Zhang       for (k=0; nextp<pnz; k++) {
65029825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
65129825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
65229825010SHong Zhang         }
653cf1a0672SHong Zhang       }
654cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
655cf1a0672SHong Zhang     }
6561634b032SHong Zhang 
657cf1a0672SHong Zhang     /* set values in C */
658cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
659cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6601634b032SHong Zhang 
661cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
662cf1a0672SHong Zhang     ca = coa + co->i[i];
663cf1a0672SHong Zhang     k  = 0;
664cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
665cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
66629825010SHong Zhang       ca[k0]        = apa_sparse[k];
66729825010SHong Zhang       apa_sparse[k] = 0.0;
668cf1a0672SHong Zhang       k++;
669cf1a0672SHong Zhang     }
6701634b032SHong Zhang 
671cf1a0672SHong Zhang     /* diagonal part of C */
672cf1a0672SHong Zhang     ca = cda + cd->i[i];
673cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
67429825010SHong Zhang       ca[k1]        = apa_sparse[k];
67529825010SHong Zhang       apa_sparse[k] = 0.0;
676cf1a0672SHong Zhang       k++;
677cf1a0672SHong Zhang     }
678cf1a0672SHong Zhang 
679cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
680cf1a0672SHong Zhang     ca = coa + co->i[i];
681cf1a0672SHong Zhang     for (; k0<conz; k0++) {
68229825010SHong Zhang       ca[k0]        = apa_sparse[k];
68329825010SHong Zhang       apa_sparse[k] = 0.0;
684cf1a0672SHong Zhang       k++;
685cf1a0672SHong Zhang     }
686cf1a0672SHong Zhang   }
687cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
688cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
689cf1a0672SHong Zhang   PetscFunctionReturn(0);
690cf1a0672SHong Zhang }
691cf1a0672SHong Zhang 
6920fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
693cf1a0672SHong Zhang #undef __FUNCT__
694b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
695b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
696cf1a0672SHong Zhang {
697cf1a0672SHong Zhang   PetscErrorCode     ierr;
698ce94432eSBarry Smith   MPI_Comm           comm;
699a7c7454dSHong Zhang   PetscMPIInt        size;
700cf1a0672SHong Zhang   Mat                Cmpi;
701cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
7020298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
703cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
704cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
705cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
706cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
707f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
708eca6b86aSHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n,Crmax;
709cf1a0672SHong Zhang   PetscReal          afill;
710cf1a0672SHong Zhang   PetscScalar        *apa;
711eca6b86aSHong Zhang   PetscTable         ta;
712cf1a0672SHong Zhang 
713cf1a0672SHong Zhang   PetscFunctionBegin;
714ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
715a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
716a7c7454dSHong Zhang 
717cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
718b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
719cf1a0672SHong Zhang 
720cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
721b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
72219f0e57aSHong Zhang 
723cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
724cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
725cf1a0672SHong Zhang 
726cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
727cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
728a7c7454dSHong Zhang   if (size > 1) {
729a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
73020e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
731968382fdSHong Zhang   } else {
73220e3dc75SHong Zhang     p_oth  = NULL;
73320e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
734a7c7454dSHong Zhang   }
735cf1a0672SHong Zhang 
736cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
737cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
738854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
739cf1a0672SHong Zhang   ptap->api = api;
740cf1a0672SHong Zhang   api[0]    = 0;
741cf1a0672SHong Zhang 
742cf1a0672SHong Zhang   /* create and initialize a linked list */
743e717af80SHong Zhang   Crmax = 6*(p_loc->rmax + (PetscInt)(1.e-2*pN)); /* expected Crmax */
744eca6b86aSHong Zhang   if (Crmax > pN) Crmax = pN;
745eca6b86aSHong Zhang   ierr = PetscTableCreate(Crmax,pN,&ta);CHKERRQ(ierr);
746eca6b86aSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
747eca6b86aSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
748eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
749eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
750eca6b86aSHong Zhang 
751eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
752cf1a0672SHong Zhang 
753cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
754cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7552205254eSKarl Rupp 
756cf1a0672SHong Zhang   current_space = free_space;
757cf1a0672SHong Zhang 
758cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
759cf1a0672SHong Zhang   for (i=0; i<am; i++) {
760cf1a0672SHong Zhang     /* diagonal portion of A */
761cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
762cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
763cf1a0672SHong Zhang       row  = *adj++;
764cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
765cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
766cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
767f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
768cf1a0672SHong Zhang     }
769cf1a0672SHong Zhang     /* off-diagonal portion of A */
770cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
771cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
772cf1a0672SHong Zhang       row  = *aoj++;
773cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
774cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
775f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
776cf1a0672SHong Zhang     }
777cf1a0672SHong Zhang 
778f84c536bSHong Zhang     apnz     = *lnk;
779cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
780cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
781cf1a0672SHong Zhang 
782cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
783cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
784cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
785cf1a0672SHong Zhang       nspacedouble++;
786cf1a0672SHong Zhang     }
787cf1a0672SHong Zhang 
788cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
789f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
790cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7912205254eSKarl Rupp 
792cf1a0672SHong Zhang     current_space->array           += apnz;
793cf1a0672SHong Zhang     current_space->local_used      += apnz;
794cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
795cf1a0672SHong Zhang   }
796cf1a0672SHong Zhang 
797cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
798cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
799854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
800cf1a0672SHong Zhang   apj  = ptap->apj;
801cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
802f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
803cf1a0672SHong Zhang 
804cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
805cf1a0672SHong Zhang   /*----------------------------------------------------*/
806cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
807cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
80833d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
809cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
810cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
811cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
812cf1a0672SHong Zhang 
81329825010SHong Zhang   /* malloc apa for assembly Cmpi */
8141795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8152205254eSKarl Rupp 
81629825010SHong Zhang   ptap->apa = apa;
817cf1a0672SHong Zhang   for (i=0; i<am; i++) {
818cf1a0672SHong Zhang     row  = i + rstart;
819cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
820cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
821cf1a0672SHong Zhang     apj += apnz;
822cf1a0672SHong Zhang   }
823cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
824cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
825cf1a0672SHong Zhang 
826cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
827cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
828f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
829cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
830cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
831cf1a0672SHong Zhang 
832cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
833cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
834cf1a0672SHong Zhang   c->ptap = ptap;
835cf1a0672SHong Zhang 
836cf1a0672SHong Zhang   *C = Cmpi;
837cf1a0672SHong Zhang 
838cf1a0672SHong Zhang   /* set MatInfo */
839118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
840cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
841cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
842cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
843cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
844cf1a0672SHong Zhang 
845cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
846cf1a0672SHong Zhang   if (api[am]) {
84757622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
84857622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
849cf1a0672SHong Zhang   } else {
850cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
851cf1a0672SHong Zhang   }
852cf1a0672SHong Zhang #endif
8531634b032SHong Zhang   PetscFunctionReturn(0);
8541634b032SHong Zhang }
855f96d37fcSHong Zhang 
856f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
857f96d37fcSHong Zhang #undef __FUNCT__
858f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
859c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
860f96d37fcSHong Zhang {
861f96d37fcSHong Zhang   PetscErrorCode ierr;
862c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
863c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
864f96d37fcSHong Zhang 
865f96d37fcSHong Zhang   PetscFunctionBegin;
866c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
867c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
868c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
869c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
870c216dbf3SHong Zhang 
871c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
872c216dbf3SHong Zhang     switch (alg) {
873c216dbf3SHong Zhang     case 1:
8742bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
875c216dbf3SHong Zhang       break;
876c216dbf3SHong Zhang     case 2:
877275476c6SMatthew G. Knepley     {
878ab1b013aSHong Zhang       Mat         Pt;
879ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
880ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
881ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
882ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
883ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
884ab1b013aSHong Zhang       ptap     = c->ptap;
885ab1b013aSHong Zhang       ptap->Pt = Pt;
886c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
887c216dbf3SHong Zhang       PetscFunctionReturn(0);
888275476c6SMatthew G. Knepley     }
889c216dbf3SHong Zhang       break;
890c216dbf3SHong Zhang     default:
8916da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
892c216dbf3SHong Zhang       break;
893c216dbf3SHong Zhang     }
894c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
895c216dbf3SHong Zhang   }
896c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
897c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
898c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
899ab1b013aSHong Zhang   PetscFunctionReturn(0);
900ab1b013aSHong Zhang }
901ab1b013aSHong Zhang 
902c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
903c216dbf3SHong Zhang #undef __FUNCT__
904c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
905c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
906c216dbf3SHong Zhang {
907c216dbf3SHong Zhang   PetscErrorCode ierr;
9082bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9092bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9102bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
911c216dbf3SHong Zhang 
912c216dbf3SHong Zhang   PetscFunctionBegin;
913c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
914c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
915f96d37fcSHong Zhang   PetscFunctionReturn(0);
916f96d37fcSHong Zhang }
917f96d37fcSHong Zhang 
9186da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
919f96d37fcSHong Zhang #undef __FUNCT__
9206da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9216da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
922f96d37fcSHong Zhang {
923c5af039cSHong Zhang   PetscErrorCode      ierr;
924c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
925497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
926c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
927c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
928d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
929497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
930e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
931c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
932ce94432eSBarry Smith   MPI_Comm            comm;
933c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
934c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
935c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
936c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
937c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
938c5af039cSHong Zhang   MPI_Status          *status;
939c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
940d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
941a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
942c5af039cSHong Zhang   Mat                 A_loc;
943c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
944f96d37fcSHong Zhang 
945f96d37fcSHong Zhang   PetscFunctionBegin;
946ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
947c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
948c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
949c5af039cSHong Zhang 
950c5af039cSHong Zhang   ptap  = c->ptap;
951c5af039cSHong Zhang   merge = ptap->merge;
952c5af039cSHong Zhang 
953c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
954c5af039cSHong Zhang   /*--------------------------------------------------------------*/
955c5af039cSHong Zhang   /* get data from symbolic products */
956c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
957854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
958c5af039cSHong Zhang 
959c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
960c5af039cSHong Zhang   owners = merge->rowmap->range;
961854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
962c5af039cSHong Zhang 
963c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
964c5af039cSHong Zhang   A_loc = ptap->A_loc;
965c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
966c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
967d6ab1dc0SHong Zhang   ai    = a_loc->i;
968d6ab1dc0SHong Zhang   aj    = a_loc->j;
969c5af039cSHong Zhang 
970854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
971c5af039cSHong Zhang 
972c5af039cSHong Zhang   for (i=0; i<am; i++) {
973e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
974d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
975d6ab1dc0SHong Zhang     adj = aj + ai[i];
976d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
977c5af039cSHong Zhang     for (j=0; j<anz; j++) {
978e745005bSHong Zhang       aval[adj[j]] = ada[j];
979c5af039cSHong Zhang     }
980c5af039cSHong Zhang 
981c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
982c5af039cSHong Zhang     /*--------------------------------------------------------------*/
983c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
984c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
985c5af039cSHong Zhang     poJ = po->j + po->i[i];
986c5af039cSHong Zhang     pA  = po->a + po->i[i];
987c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
988c5af039cSHong Zhang       row = poJ[j];
989c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
990c5af039cSHong Zhang       cj  = coj + coi[row];
991c5af039cSHong Zhang       ca  = coa + coi[row];
992c5af039cSHong Zhang       /* perform dense axpy */
993c5af039cSHong Zhang       valtmp = pA[j];
994c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
995e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
996c5af039cSHong Zhang       }
997c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
998c5af039cSHong Zhang     }
999c5af039cSHong Zhang 
1000c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1001c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1002c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1003c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1004c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1005c5af039cSHong Zhang       row = pdJ[j];
1006c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1007c5af039cSHong Zhang       cj  = bj + bi[row];
1008c5af039cSHong Zhang       ca  = ba + bi[row];
1009c5af039cSHong Zhang       /* perform dense axpy */
1010c5af039cSHong Zhang       valtmp = pA[j];
1011c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1012e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1013c5af039cSHong Zhang       }
1014c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1015c5af039cSHong Zhang     }
1016c5af039cSHong Zhang 
1017d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1018d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1019e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1020c5af039cSHong Zhang   }
1021c5af039cSHong Zhang 
1022c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1023c5af039cSHong Zhang   /*------------------------------------*/
1024c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1025c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1026c5af039cSHong Zhang   len_s  = merge->len_s;
1027c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1028c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1029c5af039cSHong Zhang 
1030dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1031c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1032c5af039cSHong Zhang     if (!len_s[proc]) continue;
1033c5af039cSHong Zhang     i    = merge->owners_co[proc];
1034c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1035c5af039cSHong Zhang     k++;
1036c5af039cSHong Zhang   }
1037c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1038c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1039c5af039cSHong Zhang 
1040c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1041c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1042c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1043c5af039cSHong Zhang 
1044c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1045c5af039cSHong Zhang   /*----------------------------------------------------*/
1046dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1047c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1048c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1049c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1050c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1051c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1052c5af039cSHong Zhang   }
1053c5af039cSHong Zhang 
1054c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1055c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1056c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1057c5af039cSHong Zhang     ba_i = ba + bi[i];
1058c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1059c5af039cSHong Zhang     /* add received vals into ba */
1060c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1061c5af039cSHong Zhang       /* i-th row */
1062c5af039cSHong Zhang       if (i == *nextrow[k]) {
1063c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1064c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1065c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1066c5af039cSHong Zhang         nextcj = 0;
1067c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1068c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1069c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1070c5af039cSHong Zhang           }
1071c5af039cSHong Zhang         }
1072c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1073c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1074c5af039cSHong Zhang       }
1075c5af039cSHong Zhang     }
1076c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1077c5af039cSHong Zhang   }
1078c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1079c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1080c5af039cSHong Zhang 
1081c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1082c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1083c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1084c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1085e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1086f96d37fcSHong Zhang   PetscFunctionReturn(0);
1087f96d37fcSHong Zhang }
1088f96d37fcSHong Zhang 
1089c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1090c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1091f96d37fcSHong Zhang #undef __FUNCT__
10922bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
10932bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1094f96d37fcSHong Zhang {
1095f96d37fcSHong Zhang   PetscErrorCode      ierr;
10964e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1097f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10980298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1099f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1100f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1101f96d37fcSHong Zhang   PetscInt            nnz;
11024e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1103497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1104f96d37fcSHong Zhang   PetscBT             lnkbt;
1105ce94432eSBarry Smith   MPI_Comm            comm;
1106f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1107f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1108f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1109f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1110f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1111f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1112f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1113f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1114d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1115f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1116438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1117f96d37fcSHong Zhang   PetscScalar         *vals;
11184e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1119f96d37fcSHong Zhang 
1120f96d37fcSHong Zhang   PetscFunctionBegin;
1121ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1122c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1123*6c4ed002SBarry 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);
1124c5af039cSHong Zhang 
1125f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1126f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1127f96d37fcSHong Zhang 
1128f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1129b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1130f96d37fcSHong Zhang 
1131f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1132f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11332205254eSKarl Rupp 
1134c5af039cSHong Zhang   ptap->A_loc = A_loc;
11352205254eSKarl Rupp 
11361c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1137d6ab1dc0SHong Zhang   ai    = a_loc->i;
1138d6ab1dc0SHong Zhang   aj    = a_loc->j;
1139f96d37fcSHong Zhang 
1140f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1141f96d37fcSHong Zhang   /*----------------------------------------------------*/
11424e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11434e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11444e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11454e938277SHong Zhang 
11464e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11474e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11484e938277SHong Zhang   poti = pot->i; potj = pot->j;
1149f96d37fcSHong Zhang 
1150f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11512205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1152854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1153f96d37fcSHong Zhang   coi[0] = 0;
1154f96d37fcSHong Zhang 
1155f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1156d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1157a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1158f96d37fcSHong Zhang   current_space = free_space;
115919f0e57aSHong Zhang 
1160f96d37fcSHong Zhang   /* create and initialize a linked list */
1161438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1162f96d37fcSHong Zhang 
1163f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1164f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1165f96d37fcSHong Zhang     ptJ = potj + poti[i];
1166f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1167f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1168d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1169d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1170f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1171d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1172f96d37fcSHong Zhang     }
11734e938277SHong Zhang     nnz = lnk[0];
1174f96d37fcSHong Zhang 
1175f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1176f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1177f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1178f96d37fcSHong Zhang       nspacedouble++;
1179f96d37fcSHong Zhang     }
1180f96d37fcSHong Zhang 
1181f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11824e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11832205254eSKarl Rupp 
1184f96d37fcSHong Zhang     current_space->array           += nnz;
1185f96d37fcSHong Zhang     current_space->local_used      += nnz;
1186f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11872205254eSKarl Rupp 
1188f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1189f96d37fcSHong Zhang   }
1190f96d37fcSHong Zhang 
1191854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1192f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11932205254eSKarl Rupp 
1194118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1195f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1196f96d37fcSHong Zhang 
1197f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1198f96d37fcSHong Zhang   /*----------------------------------------------*/
1199f96d37fcSHong Zhang   /* determine row ownership */
1200b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1201f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12022205254eSKarl Rupp 
1203f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1204f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12052205254eSKarl Rupp 
1206f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1207f96d37fcSHong Zhang   owners = merge->rowmap->range;
1208f96d37fcSHong Zhang 
1209f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12101795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1211785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12122205254eSKarl Rupp 
1213f96d37fcSHong Zhang   len_s        = merge->len_s;
1214f96d37fcSHong Zhang   merge->nsend = 0;
1215f96d37fcSHong Zhang 
1216854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1217f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1218f96d37fcSHong Zhang 
1219f96d37fcSHong Zhang   proc = 0;
1220f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1221f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1222f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1223f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1224f96d37fcSHong Zhang   }
1225f96d37fcSHong Zhang 
1226f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1227f96d37fcSHong Zhang   owners_co[0] = 0;
1228f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1229f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1230f96d37fcSHong Zhang     if (len_si[proc]) {
1231f96d37fcSHong Zhang       merge->nsend++;
1232f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1233f96d37fcSHong Zhang       len         += len_si[proc];
1234f96d37fcSHong Zhang     }
1235f96d37fcSHong Zhang   }
1236f96d37fcSHong Zhang 
1237f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12380298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1239f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1240f96d37fcSHong Zhang 
1241f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1242f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1243f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1244854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1245f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1246f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1247f96d37fcSHong Zhang     i    = owners_co[proc];
1248f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1249f96d37fcSHong Zhang     k++;
1250f96d37fcSHong Zhang   }
1251f96d37fcSHong Zhang 
1252f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1253785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1254f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1255c280ed6aSJed Brown     PetscMPIInt icompleted;
1256c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1257f96d37fcSHong Zhang   }
1258f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1259f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1260f96d37fcSHong Zhang 
1261f96d37fcSHong Zhang   /* send and recv coi */
1262f96d37fcSHong Zhang   /*-------------------*/
1263f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1264f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1265854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1266f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1267f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1268f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1269f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1270f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1271f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1272f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1273f96d37fcSHong Zhang     */
1274f96d37fcSHong Zhang     /*-------------------------------------------*/
1275f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1276f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1277f96d37fcSHong Zhang     buf_si[0]   = nrows;
1278f96d37fcSHong Zhang     buf_si_i[0] = 0;
1279f96d37fcSHong Zhang     nrows       = 0;
1280f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1281f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1282f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1283f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1284f96d37fcSHong Zhang       nrows++;
1285f96d37fcSHong Zhang     }
1286f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1287f96d37fcSHong Zhang     k++;
1288f96d37fcSHong Zhang     buf_si += len_si[proc];
1289f96d37fcSHong Zhang   }
1290f96d37fcSHong Zhang   i = merge->nrecv;
1291f96d37fcSHong Zhang   while (i--) {
1292c280ed6aSJed Brown     PetscMPIInt icompleted;
1293c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1294f96d37fcSHong Zhang   }
1295f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1296f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1297f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1298f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1299f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1300f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1301f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1302f96d37fcSHong Zhang 
1303f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1304f96d37fcSHong Zhang   /*------------------------------------------*/
1305f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1306854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1307f96d37fcSHong Zhang   bi[0] = 0;
1308f96d37fcSHong Zhang 
1309c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1310d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1311a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1312f96d37fcSHong Zhang   current_space = free_space;
1313f96d37fcSHong Zhang 
1314dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1315f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1316f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1317f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1318f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1319f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1320f96d37fcSHong Zhang   }
1321f96d37fcSHong Zhang 
13221c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1323f96d37fcSHong Zhang   rmax = 0;
1324f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1325f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1326f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1327f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1328f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1329f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1330d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1331d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1332f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1333d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1334f96d37fcSHong Zhang     }
13354e938277SHong Zhang 
1336f96d37fcSHong Zhang     /* add received col data into lnk */
1337f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1338f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1339f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1340f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13414e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1342f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1343f96d37fcSHong Zhang       }
1344f96d37fcSHong Zhang     }
13454e938277SHong Zhang     nnz = lnk[0];
1346f96d37fcSHong Zhang 
1347f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1348f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1349f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1350f96d37fcSHong Zhang       nspacedouble++;
1351f96d37fcSHong Zhang     }
1352f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13534e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1354f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13552205254eSKarl Rupp 
1356f96d37fcSHong Zhang     current_space->array           += nnz;
1357f96d37fcSHong Zhang     current_space->local_used      += nnz;
1358f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13592205254eSKarl Rupp 
1360f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1361f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1362f96d37fcSHong Zhang   }
1363f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1364f96d37fcSHong Zhang 
1365854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1366f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13672205254eSKarl Rupp 
1368118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1369f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1370d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13714e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13724e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1373f96d37fcSHong Zhang 
13741c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13751c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13761795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1377f96d37fcSHong Zhang 
1378f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1379f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
138033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1381f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1382f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1383f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1384f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1385f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1386f96d37fcSHong Zhang     row  = i + rstart;
1387f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1388f96d37fcSHong Zhang     Jptr = bj + bi[i];
1389f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1390f96d37fcSHong Zhang   }
1391f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1392f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1393f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1394f96d37fcSHong Zhang 
1395f96d37fcSHong Zhang   merge->bi        = bi;
1396f96d37fcSHong Zhang   merge->bj        = bj;
1397f96d37fcSHong Zhang   merge->coi       = coi;
1398f96d37fcSHong Zhang   merge->coj       = coj;
1399f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1400f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1401f96d37fcSHong Zhang   merge->owners_co = owners_co;
1402f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1403f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1404f96d37fcSHong Zhang 
14056da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1406c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1407c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1408f96d37fcSHong Zhang 
1409f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1410f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1411f96d37fcSHong Zhang   c->ptap     = ptap;
14120298fd71SBarry Smith   ptap->api   = NULL;
14130298fd71SBarry Smith   ptap->apj   = NULL;
1414f96d37fcSHong Zhang   ptap->merge = merge;
1415d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1416d6ab1dc0SHong Zhang 
1417d6ab1dc0SHong Zhang   *C = Cmpi;
1418d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1419d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
142057622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
142157622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1422d6ab1dc0SHong Zhang   } else {
1423d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1424d6ab1dc0SHong Zhang   }
1425d6ab1dc0SHong Zhang #endif
1426d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1427d6ab1dc0SHong Zhang }
1428d6ab1dc0SHong Zhang 
1429d6ab1dc0SHong Zhang #undef __FUNCT__
14306da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14316da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1432d6ab1dc0SHong Zhang {
1433d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1434d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1435d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1436d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1437d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1438e745005bSHong Zhang   PetscInt            *adj;
1439e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1440e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1441d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1442ce94432eSBarry Smith   MPI_Comm            comm;
1443d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1444d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1445d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1446d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1447d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1448d6ab1dc0SHong Zhang   MPI_Status          *status;
1449d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1450d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1451a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1452d6ab1dc0SHong Zhang   Mat                 A_loc;
1453d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1454d6ab1dc0SHong Zhang 
1455d6ab1dc0SHong Zhang   PetscFunctionBegin;
1456ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1457d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1458d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1459d6ab1dc0SHong Zhang 
1460d6ab1dc0SHong Zhang   ptap  = c->ptap;
1461d6ab1dc0SHong Zhang   merge = ptap->merge;
1462d6ab1dc0SHong Zhang 
1463e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1464e745005bSHong Zhang   /*------------------------------------------*/
1465d6ab1dc0SHong Zhang   /* get data from symbolic products */
1466d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1467854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1468d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1469d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14701795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1471d6ab1dc0SHong Zhang 
1472d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1473d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1474d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1475d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1476d6ab1dc0SHong Zhang   ai    = a_loc->i;
1477d6ab1dc0SHong Zhang   aj    = a_loc->j;
1478d6ab1dc0SHong Zhang 
1479d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1480d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1481d6ab1dc0SHong Zhang     adj = aj + ai[i];
1482d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1483d6ab1dc0SHong Zhang 
1484d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1485e745005bSHong Zhang     /*-------------------------------------------------------------*/
1486d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1487d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1488d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1489d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1490d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1491d6ab1dc0SHong Zhang       row = poJ[j];
1492d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1493d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1494e745005bSHong Zhang       /* perform sparse axpy */
1495e745005bSHong Zhang       nexta  = 0;
1496d6ab1dc0SHong Zhang       valtmp = pA[j];
1497e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1498e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1499e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1500e745005bSHong Zhang           nexta++;
1501d6ab1dc0SHong Zhang         }
1502e745005bSHong Zhang       }
1503e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1504d6ab1dc0SHong Zhang     }
1505d6ab1dc0SHong Zhang 
1506d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1507d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1508d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1509d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1510d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1511d6ab1dc0SHong Zhang       row = pdJ[j];
1512d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1513d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1514e745005bSHong Zhang       /* perform sparse axpy */
1515e745005bSHong Zhang       nexta  = 0;
1516d6ab1dc0SHong Zhang       valtmp = pA[j];
1517e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1518e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1519e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1520e745005bSHong Zhang           nexta++;
1521d6ab1dc0SHong Zhang         }
1522e745005bSHong Zhang       }
1523e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1524d6ab1dc0SHong Zhang     }
1525d6ab1dc0SHong Zhang   }
1526d6ab1dc0SHong Zhang 
1527d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1528d6ab1dc0SHong Zhang   /*------------------------------------*/
1529d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1530d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1531d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1532d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1533d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1534d6ab1dc0SHong Zhang 
1535dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1536d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1537d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1538d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1539d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1540d6ab1dc0SHong Zhang     k++;
1541d6ab1dc0SHong Zhang   }
1542d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1543d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1544d6ab1dc0SHong Zhang 
1545d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1546d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1547d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1548d6ab1dc0SHong Zhang 
1549d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1550d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1551dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1552d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1553e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1554d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1555d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1556d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1557d6ab1dc0SHong Zhang   }
1558d6ab1dc0SHong Zhang 
1559d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1560d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1561d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1562d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1563d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1564d6ab1dc0SHong Zhang     /* add received vals into ba */
1565d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1566d6ab1dc0SHong Zhang       /* i-th row */
1567d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1568d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1569d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1570d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1571d6ab1dc0SHong Zhang         nextcj = 0;
1572d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1573d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1574d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1575d6ab1dc0SHong Zhang           }
1576d6ab1dc0SHong Zhang         }
1577d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1578d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1579d6ab1dc0SHong Zhang       }
1580d6ab1dc0SHong Zhang     }
1581d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1582d6ab1dc0SHong Zhang   }
1583d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1584d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1585d6ab1dc0SHong Zhang 
1586d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1587d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1588d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1589d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1590d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1591d6ab1dc0SHong Zhang }
1592d6ab1dc0SHong Zhang 
1593c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
15942bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15952bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1596d6ab1dc0SHong Zhang #undef __FUNCT__
15976da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
15986da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1599d6ab1dc0SHong Zhang {
1600d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1601d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1602d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
16030298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1604d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1605d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1606d6ab1dc0SHong Zhang   PetscInt            nnz;
1607d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1608d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1609ce94432eSBarry Smith   MPI_Comm            comm;
1610d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1611d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1612d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1613d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1614d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1615d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1616d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1617d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1618d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1619d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
16204b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1621d6ab1dc0SHong Zhang   PetscScalar         *vals;
1622d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
16230acc65b4SHong Zhang   PetscTable          ta;
1624d6ab1dc0SHong Zhang 
1625d6ab1dc0SHong Zhang   PetscFunctionBegin;
1626ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1627d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1628*6c4ed002SBarry 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);
1629d6ab1dc0SHong Zhang 
1630d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1631d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1632d6ab1dc0SHong Zhang 
1633d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1634b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1635d6ab1dc0SHong Zhang 
1636d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1637d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16382205254eSKarl Rupp 
1639d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1640d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1641d6ab1dc0SHong Zhang   ai          = a_loc->i;
1642d6ab1dc0SHong Zhang   aj          = a_loc->j;
1643d6ab1dc0SHong Zhang 
1644d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1645d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1646d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1647d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1648d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1649d6ab1dc0SHong Zhang 
1650d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1651d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1652d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1653d6ab1dc0SHong Zhang 
1654d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1655d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1656d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1657854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1658d6ab1dc0SHong Zhang   coi[0] = 0;
1659d6ab1dc0SHong Zhang 
1660d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1661d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1662a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1663d6ab1dc0SHong Zhang   current_space = free_space;
166419f0e57aSHong Zhang 
1665d6ab1dc0SHong Zhang   /* create and initialize a linked list */
16660acc65b4SHong Zhang   ierr = PetscTableCreate(2*a_loc->rmax,aN,&ta);CHKERRQ(ierr);
16674b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
16680acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
16690acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1670d6ab1dc0SHong Zhang 
1671d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1672d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1673d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1674d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1675d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1676d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1677d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1678d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1679d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1680d6ab1dc0SHong Zhang     }
1681d6ab1dc0SHong Zhang     nnz = lnk[0];
1682d6ab1dc0SHong Zhang 
1683d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1684d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1685d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1686d6ab1dc0SHong Zhang       nspacedouble++;
1687d6ab1dc0SHong Zhang     }
1688d6ab1dc0SHong Zhang 
1689d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1690d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16912205254eSKarl Rupp 
1692d6ab1dc0SHong Zhang     current_space->array           += nnz;
1693d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1694d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16952205254eSKarl Rupp 
1696d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1697d6ab1dc0SHong Zhang   }
1698d6ab1dc0SHong Zhang 
1699854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1700d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17010acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
17022205254eSKarl Rupp 
1703118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1704d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1705d6ab1dc0SHong Zhang 
1706d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1707d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1708d6ab1dc0SHong Zhang   /* determine row ownership */
1709b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1710d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17112205254eSKarl Rupp 
1712d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1713d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17142205254eSKarl Rupp 
1715d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1716d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1717d6ab1dc0SHong Zhang 
1718d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17191795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1720785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17212205254eSKarl Rupp 
1722d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1723d6ab1dc0SHong Zhang   merge->nsend = 0;
1724d6ab1dc0SHong Zhang 
1725854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1726d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1727d6ab1dc0SHong Zhang 
1728d6ab1dc0SHong Zhang   proc = 0;
1729d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1730d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1731d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1732d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1733d6ab1dc0SHong Zhang   }
1734d6ab1dc0SHong Zhang 
1735d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1736d6ab1dc0SHong Zhang   owners_co[0] = 0;
1737d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1738d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1739d6ab1dc0SHong Zhang     if (len_si[proc]) {
1740d6ab1dc0SHong Zhang       merge->nsend++;
1741d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1742d6ab1dc0SHong Zhang       len         += len_si[proc];
1743d6ab1dc0SHong Zhang     }
1744d6ab1dc0SHong Zhang   }
1745d6ab1dc0SHong Zhang 
1746d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17470298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1748d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1749d6ab1dc0SHong Zhang 
1750d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1751d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1752d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1753854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1754d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1755d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1756d6ab1dc0SHong Zhang     i    = owners_co[proc];
1757d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1758d6ab1dc0SHong Zhang     k++;
1759d6ab1dc0SHong Zhang   }
1760d6ab1dc0SHong Zhang 
1761d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1762785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1764c280ed6aSJed Brown     PetscMPIInt icompleted;
1765c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1766d6ab1dc0SHong Zhang   }
1767d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1768d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1769d6ab1dc0SHong Zhang 
17700acc65b4SHong Zhang   /* add received column indices into table to update Armax */
17710acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
17720acc65b4SHong Zhang     Jptr = buf_rj[k];
17730acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
17740acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
17750acc65b4SHong Zhang     }
17760acc65b4SHong Zhang   }
17770acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
17780acc65b4SHong Zhang 
1779d6ab1dc0SHong Zhang   /* send and recv coi */
1780d6ab1dc0SHong Zhang   /*-------------------*/
1781d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1782d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1783854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1784d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1785d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1786d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1787d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1788d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1789d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1790d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1791d6ab1dc0SHong Zhang     */
1792d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1793d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1794d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1795d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1796d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1797d6ab1dc0SHong Zhang     nrows       = 0;
1798d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1799d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1800d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1801d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1802d6ab1dc0SHong Zhang       nrows++;
1803d6ab1dc0SHong Zhang     }
1804d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1805d6ab1dc0SHong Zhang     k++;
1806d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1807d6ab1dc0SHong Zhang   }
1808d6ab1dc0SHong Zhang   i = merge->nrecv;
1809d6ab1dc0SHong Zhang   while (i--) {
1810c280ed6aSJed Brown     PetscMPIInt icompleted;
1811c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1812d6ab1dc0SHong Zhang   }
1813d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1815d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1816d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1817d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1818d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1820d6ab1dc0SHong Zhang 
1821d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1822d6ab1dc0SHong Zhang   /*------------------------------------------*/
1823d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1824854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1825d6ab1dc0SHong Zhang   bi[0] = 0;
1826d6ab1dc0SHong Zhang 
1827d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1828d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1829a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1830d6ab1dc0SHong Zhang   current_space = free_space;
1831d6ab1dc0SHong Zhang 
1832dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1833d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1834d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1835d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1836d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18372205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1838d6ab1dc0SHong Zhang   }
1839d6ab1dc0SHong Zhang 
18400acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1841d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1842d6ab1dc0SHong Zhang   rmax = 0;
1843d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1844d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1845d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1846d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1847d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1848d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1849d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1850d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1851d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1852d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1853d6ab1dc0SHong Zhang     }
1854d6ab1dc0SHong Zhang 
1855d6ab1dc0SHong Zhang     /* add received col data into lnk */
1856d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1857d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1858d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1859d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1860d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1861d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1862d6ab1dc0SHong Zhang       }
1863d6ab1dc0SHong Zhang     }
1864d6ab1dc0SHong Zhang     nnz = lnk[0];
1865d6ab1dc0SHong Zhang 
1866d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1867d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1868d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1869d6ab1dc0SHong Zhang       nspacedouble++;
1870d6ab1dc0SHong Zhang     }
1871d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1872d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1873d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18742205254eSKarl Rupp 
1875d6ab1dc0SHong Zhang     current_space->array           += nnz;
1876d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1877d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18782205254eSKarl Rupp 
1879d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1880d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1881d6ab1dc0SHong Zhang   }
1882f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1883d6ab1dc0SHong Zhang 
1884854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1885d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1886118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1887d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1888d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
18890acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
18900acc65b4SHong Zhang 
1891d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1892d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1893d6ab1dc0SHong Zhang 
1894d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1895d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18961795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1897d6ab1dc0SHong Zhang 
1898d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1899d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
190033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1901d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1902d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1903d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1904d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1905d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1906d6ab1dc0SHong Zhang     row  = i + rstart;
1907d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1908d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1909d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1910d6ab1dc0SHong Zhang   }
1911d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1912d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1913d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1914d6ab1dc0SHong Zhang 
1915d6ab1dc0SHong Zhang   merge->bi        = bi;
1916d6ab1dc0SHong Zhang   merge->bj        = bj;
1917d6ab1dc0SHong Zhang   merge->coi       = coi;
1918d6ab1dc0SHong Zhang   merge->coj       = coj;
1919d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1920d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1921d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1922d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1923d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1924d6ab1dc0SHong Zhang 
19256da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1926d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1927c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1928d6ab1dc0SHong Zhang 
1929d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1930d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19312205254eSKarl Rupp 
1932d6ab1dc0SHong Zhang   c->ptap     = ptap;
19330298fd71SBarry Smith   ptap->api   = NULL;
19340298fd71SBarry Smith   ptap->apj   = NULL;
1935d6ab1dc0SHong Zhang   ptap->merge = merge;
1936d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19370298fd71SBarry Smith   ptap->apa   = NULL;
1938f96d37fcSHong Zhang 
1939f96d37fcSHong Zhang   *C = Cmpi;
1940f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1941f96d37fcSHong Zhang   if (bi[pn] != 0) {
194257622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
194357622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1944f96d37fcSHong Zhang   } else {
1945f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1946f96d37fcSHong Zhang   }
1947f96d37fcSHong Zhang #endif
1948f96d37fcSHong Zhang   PetscFunctionReturn(0);
1949f96d37fcSHong Zhang }
1950