xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 904d1e7050d8278cb5f67df70e5edfa4eaab048b)
1be1d678aSKris Buschelman 
22499ec78SHong Zhang /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
42499ec78SHong Zhang           C = A * B
52499ec78SHong Zhang */
6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
8c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>
11af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
121634b032SHong Zhang 
135e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
145e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
155e5acdf2Sstefano_zampini #endif
165e5acdf2Sstefano_zampini 
17150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
182499ec78SHong Zhang {
192499ec78SHong Zhang   PetscErrorCode ierr;
205e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
21a9d7e0ccSandi selinger   const char     *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"};
22a9d7e0ccSandi selinger   PetscInt       nalg = 4;
235e5acdf2Sstefano_zampini #else
24a9d7e0ccSandi selinger   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
25a9d7e0ccSandi selinger   PetscInt       nalg = 3;
265e5acdf2Sstefano_zampini #endif
2760106aa7SHong Zhang   PetscInt       alg = 1; /* set nonscalable algorithm as default */
28e55be12dSHong Zhang   MPI_Comm       comm;
295e68f438SHong Zhang   PetscBool      flg;
302499ec78SHong Zhang 
312499ec78SHong Zhang   PetscFunctionBegin;
322499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
33e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
346c4ed002SBarry Smith     if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
35e55be12dSHong Zhang 
360fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
3768eacb73SHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
385e68f438SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr);
390fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
400fc8cf34SHong Zhang 
4160106aa7SHong Zhang     if (!flg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
425e68f438SHong Zhang       MatInfo     Ainfo,Binfo;
4360106aa7SHong Zhang       PetscInt    nz_local;
4460106aa7SHong Zhang       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
455e68f438SHong Zhang 
465e68f438SHong Zhang       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
475e68f438SHong Zhang       ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr);
485e68f438SHong Zhang       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
495e68f438SHong Zhang 
5060106aa7SHong Zhang       if (B->cmap->N > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
5160106aa7SHong Zhang       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
5260106aa7SHong Zhang 
5360106aa7SHong Zhang       if (alg_scalable) {
5460106aa7SHong Zhang         alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
5560106aa7SHong Zhang         ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,fill*nz_local);CHKERRQ(ierr);
5660106aa7SHong Zhang       }
575e68f438SHong Zhang     }
585e68f438SHong Zhang 
593ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
600fc8cf34SHong Zhang     switch (alg) {
610fc8cf34SHong Zhang     case 1:
620fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
630fc8cf34SHong Zhang       break;
645e5acdf2Sstefano_zampini     case 2:
65a9d7e0ccSandi selinger       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);CHKERRQ(ierr);
66a9d7e0ccSandi selinger       break;
67a9d7e0ccSandi selinger #if defined(PETSC_HAVE_HYPRE)
68a9d7e0ccSandi selinger     case 3:
695e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
705e5acdf2Sstefano_zampini       break;
715e5acdf2Sstefano_zampini #endif
720fc8cf34SHong Zhang     default:
73b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
740fc8cf34SHong Zhang       break;
750fc8cf34SHong Zhang     }
763ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
772499ec78SHong Zhang   }
783ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
79598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
803ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
812499ec78SHong Zhang   PetscFunctionReturn(0);
822499ec78SHong Zhang }
832499ec78SHong Zhang 
84a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
85598bc09dSHong Zhang {
86598bc09dSHong Zhang   PetscErrorCode ierr;
87598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
88598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
89598bc09dSHong Zhang 
90598bc09dSHong Zhang   PetscFunctionBegin;
91b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
92598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
93a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
94a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
95ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
96a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
97a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
98598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
99598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
100598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
101598bc09dSHong Zhang   PetscFunctionReturn(0);
102598bc09dSHong Zhang }
103598bc09dSHong Zhang 
104a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
1054ae0847bSHong Zhang {
1064ae0847bSHong Zhang   PetscErrorCode ierr;
1074ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
1084ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
1094ae0847bSHong Zhang 
1104ae0847bSHong Zhang   PetscFunctionBegin;
1114ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
11226fbe8dcSKarl Rupp 
1134ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
1144ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
1154ae0847bSHong Zhang   PetscFunctionReturn(0);
1164ae0847bSHong Zhang }
1174ae0847bSHong Zhang 
118f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
119598bc09dSHong Zhang {
120598bc09dSHong Zhang   PetscErrorCode ierr;
1214ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
122598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
123598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
124c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
125598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
1269ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
127c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
128598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
129c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
13029825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
131598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1329467f45fSHong Zhang   MPI_Comm       comm;
1339467f45fSHong Zhang   PetscMPIInt    size;
134598bc09dSHong Zhang 
135598bc09dSHong Zhang   PetscFunctionBegin;
1369467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1379467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1389467f45fSHong Zhang 
139a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
140598bc09dSHong Zhang   /*-----------------------------------------------------*/
141a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
142b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
143a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
144598bc09dSHong Zhang 
145598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
146598bc09dSHong Zhang   /*----------------------------------------------------------*/
147598bc09dSHong Zhang   /* get data from symbolic products */
148a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
149c12557a1SHong Zhang   p_oth = NULL;
1509467f45fSHong Zhang   if (size >1) {
1519467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1529467f45fSHong Zhang   }
153598bc09dSHong Zhang 
154598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
155598bc09dSHong Zhang   apa = ptap->apa;
156598bc09dSHong Zhang 
15729825010SHong Zhang   api = ptap->api;
15829825010SHong Zhang   apj = ptap->apj;
159598bc09dSHong Zhang   for (i=0; i<cm; i++) {
160c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
161e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
162598bc09dSHong Zhang 
163598bc09dSHong Zhang     /* set values in C */
164598bc09dSHong Zhang     apJ  = apj + api[i];
165598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
166598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
167598bc09dSHong Zhang 
168598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
169598bc09dSHong Zhang     ca = coa + co->i[i];
170598bc09dSHong Zhang     k  = 0;
171598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
172598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
173598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
1745cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
175598bc09dSHong Zhang     }
176598bc09dSHong Zhang 
177598bc09dSHong Zhang     /* diagonal part of C */
178598bc09dSHong Zhang     ca = cda + cd->i[i];
179598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
180598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
1815cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
182598bc09dSHong Zhang     }
183598bc09dSHong Zhang 
184598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
185598bc09dSHong Zhang     ca = coa + co->i[i];
186598bc09dSHong Zhang     for (; k0<conz; k0++) {
187598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
1885cab4f04SHong Zhang       apa[apJ[k++]] = 0.0;
189598bc09dSHong Zhang     }
190598bc09dSHong Zhang   }
191598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
193598bc09dSHong Zhang   PetscFunctionReturn(0);
194598bc09dSHong Zhang }
195598bc09dSHong Zhang 
1960fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
197598bc09dSHong Zhang {
198598bc09dSHong Zhang   PetscErrorCode     ierr;
199ce94432eSBarry Smith   MPI_Comm           comm;
2009467f45fSHong Zhang   PetscMPIInt        size;
201bfcd1a73SHong Zhang   Mat                Cmpi;
202598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
2030298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2044ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
205bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2067dbaa2c6Sandi selinger   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2074ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
208d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
2095cab4f04SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
210598bc09dSHong Zhang   PetscBT            lnkbt;
211598bc09dSHong Zhang   PetscScalar        *apa;
212bfcd1a73SHong Zhang   PetscReal          afill;
213598bc09dSHong Zhang 
214598bc09dSHong Zhang   PetscFunctionBegin;
215ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2169467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2179467f45fSHong Zhang 
218598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
219b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
220598bc09dSHong Zhang 
221598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
222b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
22319f0e57aSHong Zhang 
224598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
225a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
226598bc09dSHong Zhang 
227a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
228598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2299467f45fSHong Zhang   if (size > 1) {
2309467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
231598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
23220e3dc75SHong Zhang   } else {
23320e3dc75SHong Zhang     p_oth = NULL;
23420e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2359467f45fSHong Zhang   }
236598bc09dSHong Zhang 
237598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
238598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
239854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
240a1a4e84aSHong Zhang   ptap->api = api;
241598bc09dSHong Zhang   api[0]    = 0;
242598bc09dSHong Zhang 
2435cab4f04SHong Zhang   /* create and initialize a linked list */
2445cab4f04SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
245598bc09dSHong Zhang 
246bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
247f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
248598bc09dSHong Zhang   current_space = free_space;
249598bc09dSHong Zhang 
250598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
251598bc09dSHong Zhang   for (i=0; i<am; i++) {
252598bc09dSHong Zhang     /* diagonal portion of A */
253598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
254598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
255598bc09dSHong Zhang       row  = *adj++;
256598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
257598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
258598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2591fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
260598bc09dSHong Zhang     }
261598bc09dSHong Zhang     /* off-diagonal portion of A */
262598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
263598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
264598bc09dSHong Zhang       row  = *aoj++;
265598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
266598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2671fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
268598bc09dSHong Zhang     }
269598bc09dSHong Zhang 
270d14fa243SHong Zhang     apnz     = lnk[0];
271598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
272598bc09dSHong Zhang 
273598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
274598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
275f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
276598bc09dSHong Zhang       nspacedouble++;
277598bc09dSHong Zhang     }
278598bc09dSHong Zhang 
279598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2801fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
281598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2822205254eSKarl Rupp 
283598bc09dSHong Zhang     current_space->array           += apnz;
284598bc09dSHong Zhang     current_space->local_used      += apnz;
285598bc09dSHong Zhang     current_space->local_remaining -= apnz;
286598bc09dSHong Zhang   }
287598bc09dSHong Zhang 
288598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
289598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
290854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
291a1a4e84aSHong Zhang   apj  = ptap->apj;
292a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
293598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
294598bc09dSHong Zhang 
295f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2961795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2972205254eSKarl Rupp 
298f84c536bSHong Zhang   ptap->apa = apa;
299f84c536bSHong Zhang 
300bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
301598bc09dSHong Zhang   /*----------------------------------------------------*/
302bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
303bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
30433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
305a2f3521dSMark F. Adams 
306bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
307bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
30817f75172Sandi selinger 
309*904d1e70Sandi selinger   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr);
310bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
311bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3127dbaa2c6Sandi selinger   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
313598bc09dSHong Zhang 
314bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
315bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
316f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
317bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
318bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
319598bc09dSHong Zhang 
320bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
321bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
322598bc09dSHong Zhang   c->ptap = ptap;
323598bc09dSHong Zhang 
324bfcd1a73SHong Zhang   *C = Cmpi;
325bfcd1a73SHong Zhang 
326bfcd1a73SHong Zhang   /* set MatInfo */
327118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
328bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
329bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
330bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
331bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
332bfcd1a73SHong Zhang 
333bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
334bfcd1a73SHong Zhang   if (api[am]) {
33557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
33657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
337bfcd1a73SHong Zhang   } else {
338bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
339bfcd1a73SHong Zhang   }
340bfcd1a73SHong Zhang #endif
341598bc09dSHong Zhang   PetscFunctionReturn(0);
342598bc09dSHong Zhang }
343598bc09dSHong Zhang 
344150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3459123193aSHong Zhang {
3469123193aSHong Zhang   PetscErrorCode ierr;
3479123193aSHong Zhang 
3489123193aSHong Zhang   PetscFunctionBegin;
3499123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3503ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3519123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3523ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3539123193aSHong Zhang   }
3543ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3559123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3563ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3579123193aSHong Zhang   PetscFunctionReturn(0);
3589123193aSHong Zhang }
3599123193aSHong Zhang 
3604324174eSBarry Smith typedef struct {
3614324174eSBarry Smith   Mat         workB;
3624324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3634324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3644324174eSBarry Smith } MPIAIJ_MPIDense;
3654324174eSBarry Smith 
36696b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3674324174eSBarry Smith {
3684324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3694324174eSBarry Smith   PetscErrorCode  ierr;
3704324174eSBarry Smith 
3714324174eSBarry Smith   PetscFunctionBegin;
3726bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3734324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3744324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3754324174eSBarry Smith   PetscFunctionReturn(0);
3764324174eSBarry Smith }
3774324174eSBarry Smith 
378fd4e9aacSBarry Smith /*
379fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
380fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
381fd4e9aacSBarry Smith 
382fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
383fd4e9aacSBarry Smith */
384fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
385fd4e9aacSBarry Smith {
386fd4e9aacSBarry Smith   PetscErrorCode         ierr;
387fd4e9aacSBarry Smith   PetscBool              flg;
388fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
389fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
390fd4e9aacSBarry Smith   PetscContainer         container;
391fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
392fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
393fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
394fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
395fd4e9aacSBarry Smith 
396fd4e9aacSBarry Smith   PetscFunctionBegin;
397fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
398fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
399fd4e9aacSBarry Smith 
400fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
401fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
402fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
403fd4e9aacSBarry Smith 
404fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
405fd4e9aacSBarry Smith 
406fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
407fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
408fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
409fd4e9aacSBarry Smith   /* Create work arrays needed */
410fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
411fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
412fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
413fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
414fd4e9aacSBarry Smith 
415fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
416fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
417fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
418fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
419fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
420fd4e9aacSBarry Smith 
421fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
422fd4e9aacSBarry Smith   PetscFunctionReturn(0);
423fd4e9aacSBarry Smith }
424fd4e9aacSBarry Smith 
4259123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4269123193aSHong Zhang {
4279123193aSHong Zhang   PetscErrorCode         ierr;
428f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
429d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
430bf0cc555SLisandro Dalcin   PetscContainer         container;
4314324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4324324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4334324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4344324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
435d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4369123193aSHong Zhang 
4379123193aSHong Zhang   PetscFunctionBegin;
438ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
439d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
44033d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
441cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4420298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
443cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
444cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4452205254eSKarl Rupp 
446d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
447f32d5d43SBarry Smith 
448b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
449f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4500298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4514324174eSBarry Smith   /* Create work arrays needed */
452dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
453dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
454dcca6d9dSJed Brown                       from->n,&contents->rwaits,
455dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4564324174eSBarry Smith 
457ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
458bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
45996b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
460bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
461bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4629123193aSHong Zhang   PetscFunctionReturn(0);
4639123193aSHong Zhang }
4649123193aSHong Zhang 
465f32d5d43SBarry Smith /*
4662636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4672636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
468f32d5d43SBarry Smith */
4694324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
470f32d5d43SBarry Smith {
471f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
472f32d5d43SBarry Smith   PetscErrorCode         ierr;
473f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
474f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
475f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
476f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
477f32d5d43SBarry Smith   PetscInt               i,j,k;
478aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
479aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
480f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
481ce94432eSBarry Smith   MPI_Comm               comm;
482d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
483f32d5d43SBarry Smith   MPI_Status             status;
4844324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
485bf0cc555SLisandro Dalcin   PetscContainer         container;
4864324174eSBarry Smith   Mat                    workB;
487f32d5d43SBarry Smith 
488f32d5d43SBarry Smith   PetscFunctionBegin;
489ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
490bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
49129825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
492bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4934324174eSBarry Smith 
4944324174eSBarry Smith   workB = *outworkB = contents->workB;
495cf1a0672SHong Zhang   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
496f32d5d43SBarry Smith   sindices = to->indices;
497f32d5d43SBarry Smith   sstarts  = to->starts;
498f32d5d43SBarry Smith   sprocs   = to->procs;
4994324174eSBarry Smith   swaits   = contents->swaits;
5004324174eSBarry Smith   svalues  = contents->svalues;
501f32d5d43SBarry Smith 
502f32d5d43SBarry Smith   rindices = from->indices;
503f32d5d43SBarry Smith   rstarts  = from->starts;
504f32d5d43SBarry Smith   rprocs   = from->procs;
5054324174eSBarry Smith   rwaits   = contents->rwaits;
5064324174eSBarry Smith   rvalues  = contents->rvalues;
507f32d5d43SBarry Smith 
5088c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5098c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
510f32d5d43SBarry Smith 
511f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
512f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
513f32d5d43SBarry Smith   }
514f32d5d43SBarry Smith 
515f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
516f32d5d43SBarry Smith     /* pack a message at a time */
517f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
518f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
519f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5202636ff26SBarry Smith       }
5212636ff26SBarry Smith     }
522f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
523f32d5d43SBarry Smith   }
5242636ff26SBarry Smith 
525f32d5d43SBarry Smith   nrecvs = from->n;
526f32d5d43SBarry Smith   while (nrecvs) {
527f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
528f32d5d43SBarry Smith     nrecvs--;
529f32d5d43SBarry Smith     /* unpack a message at a time */
530f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
531f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
532f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5332636ff26SBarry Smith       }
5342636ff26SBarry Smith     }
535f32d5d43SBarry Smith   }
536cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
537f32d5d43SBarry Smith 
5388c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5398c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
540f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
541f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
542f32d5d43SBarry Smith   PetscFunctionReturn(0);
543f32d5d43SBarry Smith }
544f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
545f32d5d43SBarry Smith 
5469123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5479123193aSHong Zhang {
5489123193aSHong Zhang   PetscErrorCode ierr;
549f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
550f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
551f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
552f32d5d43SBarry Smith   Mat            workB;
5539123193aSHong Zhang 
5549123193aSHong Zhang   PetscFunctionBegin;
555f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
556f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
557f32d5d43SBarry Smith 
558f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5594324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
560f32d5d43SBarry Smith 
561f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
562f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5639123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5649123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5659123193aSHong Zhang   PetscFunctionReturn(0);
5669123193aSHong Zhang }
567cf1a0672SHong Zhang 
568f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5691634b032SHong Zhang {
5701634b032SHong Zhang   PetscErrorCode ierr;
571cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
572cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
573cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
574cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
575cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
576cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
577cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
57829825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
57929825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
580cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
58129825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
582cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
58329825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
58429825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
585a7c7454dSHong Zhang   MPI_Comm       comm;
586a7c7454dSHong Zhang   PetscMPIInt    size;
5871634b032SHong Zhang 
5881634b032SHong Zhang   PetscFunctionBegin;
589a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
590a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
591a7c7454dSHong Zhang 
592cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
593cf1a0672SHong Zhang   /*-----------------------------------------------------*/
594cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
595b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
596cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
597cf1a0672SHong Zhang 
598cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
599cf1a0672SHong Zhang   /*----------------------------------------------------------*/
600cf1a0672SHong Zhang   /* get data from symbolic products */
601cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
602cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
603a7c7454dSHong Zhang   if (size >1) {
604a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
605cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
60620e3dc75SHong Zhang   } else {
60720e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
608a7c7454dSHong Zhang   }
609cf1a0672SHong Zhang 
61029825010SHong Zhang   api = ptap->api;
61129825010SHong Zhang   apj = ptap->apj;
612cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
61329825010SHong Zhang     apJ = apj + api[i];
61429825010SHong Zhang 
615cf1a0672SHong Zhang     /* diagonal portion of A */
616cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
617cf1a0672SHong Zhang     adj = ad->j + adi[i];
618cf1a0672SHong Zhang     ada = ad->a + adi[i];
619cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
620cf1a0672SHong Zhang       row = adj[j];
621cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
622cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
623cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
62429825010SHong Zhang       /* perform sparse axpy */
625cf1a0672SHong Zhang       valtmp = ada[j];
62629825010SHong Zhang       nextp  = 0;
62729825010SHong Zhang       for (k=0; nextp<pnz; k++) {
62829825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
62929825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63029825010SHong Zhang         }
6311634b032SHong Zhang       }
632cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
633cf1a0672SHong Zhang     }
6341634b032SHong Zhang 
635cf1a0672SHong Zhang     /* off-diagonal portion of A */
636cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
637cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
638cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
639cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
640cf1a0672SHong Zhang       row = aoj[j];
641cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
642cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
643cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
64429825010SHong Zhang       /* perform sparse axpy */
645cf1a0672SHong Zhang       valtmp = aoa[j];
64629825010SHong Zhang       nextp  = 0;
64729825010SHong Zhang       for (k=0; nextp<pnz; k++) {
64829825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
64929825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
65029825010SHong Zhang         }
651cf1a0672SHong Zhang       }
652cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
653cf1a0672SHong Zhang     }
6541634b032SHong Zhang 
655cf1a0672SHong Zhang     /* set values in C */
656cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
657cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6581634b032SHong Zhang 
659cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
660cf1a0672SHong Zhang     ca = coa + co->i[i];
661cf1a0672SHong Zhang     k  = 0;
662cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
663cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
66429825010SHong Zhang       ca[k0]        = apa_sparse[k];
66529825010SHong Zhang       apa_sparse[k] = 0.0;
666cf1a0672SHong Zhang       k++;
667cf1a0672SHong Zhang     }
6681634b032SHong Zhang 
669cf1a0672SHong Zhang     /* diagonal part of C */
670cf1a0672SHong Zhang     ca = cda + cd->i[i];
671cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
67229825010SHong Zhang       ca[k1]        = apa_sparse[k];
67329825010SHong Zhang       apa_sparse[k] = 0.0;
674cf1a0672SHong Zhang       k++;
675cf1a0672SHong Zhang     }
676cf1a0672SHong Zhang 
677cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
678cf1a0672SHong Zhang     ca = coa + co->i[i];
679cf1a0672SHong Zhang     for (; k0<conz; k0++) {
68029825010SHong Zhang       ca[k0]        = apa_sparse[k];
68129825010SHong Zhang       apa_sparse[k] = 0.0;
682cf1a0672SHong Zhang       k++;
683cf1a0672SHong Zhang     }
684cf1a0672SHong Zhang   }
685cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
687cf1a0672SHong Zhang   PetscFunctionReturn(0);
688cf1a0672SHong Zhang }
689cf1a0672SHong Zhang 
6900fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
691b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
692cf1a0672SHong Zhang {
693cf1a0672SHong Zhang   PetscErrorCode     ierr;
694ce94432eSBarry Smith   MPI_Comm           comm;
695a7c7454dSHong Zhang   PetscMPIInt        size;
696cf1a0672SHong Zhang   Mat                Cmpi;
697cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6980298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
699cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
700cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
7017dbaa2c6Sandi selinger   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
702cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
7038a210135Sandi selinger   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
704e5ab8c0aSKarl Rupp   PetscInt           am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
705cf1a0672SHong Zhang   PetscReal          afill;
706cf1a0672SHong Zhang   PetscScalar        *apa;
707cf1a0672SHong Zhang 
708cf1a0672SHong Zhang   PetscFunctionBegin;
709ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
710a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
711a7c7454dSHong Zhang 
712cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
713b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
714cf1a0672SHong Zhang 
715cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
716b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
71719f0e57aSHong Zhang 
718cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
719cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
720cf1a0672SHong Zhang 
721cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
722cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
723a7c7454dSHong Zhang   if (size > 1) {
724a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
72520e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
726968382fdSHong Zhang   } else {
72720e3dc75SHong Zhang     p_oth  = NULL;
72820e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
729a7c7454dSHong Zhang   }
730cf1a0672SHong Zhang 
731cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
732cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
733854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
734cf1a0672SHong Zhang   ptap->api = api;
735cf1a0672SHong Zhang   api[0]    = 0;
736cf1a0672SHong Zhang 
7378a210135Sandi selinger   ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr);
738cf1a0672SHong Zhang 
739cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
740f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
741cf1a0672SHong Zhang   current_space = free_space;
742cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
743cf1a0672SHong Zhang   for (i=0; i<am; i++) {
744cf1a0672SHong Zhang     /* diagonal portion of A */
745cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
746cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
747cf1a0672SHong Zhang       row  = *adj++;
748cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
749cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
7508a210135Sandi selinger       /* Expand list if it is not long enough */
7518328dd82Sandi selinger       if (pnz+apnz_max > lsize) {
7528328dd82Sandi selinger         lsize = pnz+apnz_max;
75326465080Sandi selinger         ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr);
7548a210135Sandi selinger       }
755cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
756f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
7578a210135Sandi selinger       apnz     = *lnk; /* The first element in the list is the number of items in the list */
7588a210135Sandi selinger       api[i+1] = api[i] + apnz;
7598a210135Sandi selinger       if (apnz > apnz_max) apnz_max = apnz;
760cf1a0672SHong Zhang     }
761cf1a0672SHong Zhang     /* off-diagonal portion of A */
762cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
763cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
764cf1a0672SHong Zhang       row  = *aoj++;
765cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
766cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
7678a210135Sandi selinger       /* Expand list if it is not long enough */
7688328dd82Sandi selinger       if (pnz+apnz_max > lsize) {
7698328dd82Sandi selinger         lsize = pnz + apnz_max;
77026465080Sandi selinger         ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr);
7718a210135Sandi selinger       }
7728a210135Sandi selinger       /* add non-zero cols of P into the sorted linked list lnk */
773f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
7748a210135Sandi selinger       apnz     = *lnk;  /* The first element in the list is the number of items in the list */
7758a210135Sandi selinger       api[i+1] = api[i] + apnz;
7768a210135Sandi selinger       if (apnz > apnz_max) apnz_max = apnz;
777cf1a0672SHong Zhang     }
778f84c536bSHong Zhang     apnz     = *lnk;
779cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
7808a210135Sandi selinger     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) {
784f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(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 
81229825010SHong Zhang   /* malloc apa for assembly Cmpi */
8131795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
81429825010SHong Zhang   ptap->apa = apa;
81517f75172Sandi selinger 
816*904d1e70Sandi selinger   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr);
817cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
818cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8197dbaa2c6Sandi selinger   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
820cf1a0672SHong Zhang 
821cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
822cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
823f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
824cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
825cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
826cf1a0672SHong Zhang 
827cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
828cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
829cf1a0672SHong Zhang   c->ptap = ptap;
830cf1a0672SHong Zhang   *C = Cmpi;
831cf1a0672SHong Zhang 
832cf1a0672SHong Zhang   /* set MatInfo */
833118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
834cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
835cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
836cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
837cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
838cf1a0672SHong Zhang 
839cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
840cf1a0672SHong Zhang   if (api[am]) {
84157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
84257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
843cf1a0672SHong Zhang   } else {
844cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
845cf1a0672SHong Zhang   }
846cf1a0672SHong Zhang #endif
8471634b032SHong Zhang   PetscFunctionReturn(0);
8481634b032SHong Zhang }
849f96d37fcSHong Zhang 
850a9d7e0ccSandi selinger /* This function is needed for the seqMPI matrix-matrix multiplication.  */
851a9d7e0ccSandi selinger /* Three input arrays are merged to one output array. The size of the    */
852a9d7e0ccSandi selinger /* output array is also output. Duplicate entries only show up once.     */
853a9d7e0ccSandi selinger static void Merge3SortedArrays(PetscInt  size1, PetscInt *in1,
854a9d7e0ccSandi selinger                                PetscInt  size2, PetscInt *in2,
855a9d7e0ccSandi selinger                                PetscInt  size3, PetscInt *in3,
856a9d7e0ccSandi selinger                                PetscInt *size4, PetscInt *out)
857a9d7e0ccSandi selinger {
858a9d7e0ccSandi selinger   int i = 0, j = 0, k = 0, l = 0;
859a9d7e0ccSandi selinger 
860a9d7e0ccSandi selinger   /* Traverse all three arrays */
861a9d7e0ccSandi selinger   while (i<size1 && j<size2 && k<size3) {
862a9d7e0ccSandi selinger     if (in1[i] < in2[j] && in1[i] < in3[k]) {
863a9d7e0ccSandi selinger       out[l++] = in1[i++];
864a9d7e0ccSandi selinger     }
865a9d7e0ccSandi selinger     else if(in2[j] < in1[i] && in2[j] < in3[k]) {
866a9d7e0ccSandi selinger       out[l++] = in2[j++];
867a9d7e0ccSandi selinger     }
868a9d7e0ccSandi selinger     else if(in3[k] < in1[i] && in3[k] < in2[j]) {
869a9d7e0ccSandi selinger       out[l++] = in3[k++];
870a9d7e0ccSandi selinger     }
871a9d7e0ccSandi selinger     else if(in1[i] == in2[j] && in1[i] < in3[k]) {
872a9d7e0ccSandi selinger       out[l++] = in1[i];
873a9d7e0ccSandi selinger       i++, j++;
874a9d7e0ccSandi selinger     }
875a9d7e0ccSandi selinger     else if(in1[i] == in3[k] && in1[i] < in2[j]) {
876a9d7e0ccSandi selinger       out[l++] = in1[i];
877a9d7e0ccSandi selinger       i++, k++;
878a9d7e0ccSandi selinger     }
879a9d7e0ccSandi selinger     else if(in3[k] == in2[j] && in2[j] < in1[i])  {
880a9d7e0ccSandi selinger       out[l++] = in2[j];
881a9d7e0ccSandi selinger       k++, j++;
882a9d7e0ccSandi selinger     }
883a9d7e0ccSandi selinger     else if(in1[i] == in2[j] && in1[i] == in3[k]) {
884a9d7e0ccSandi selinger       out[l++] = in1[i];
885a9d7e0ccSandi selinger       i++, j++, k++;
886a9d7e0ccSandi selinger     }
887a9d7e0ccSandi selinger   }
888a9d7e0ccSandi selinger 
889a9d7e0ccSandi selinger   /* Traverse two remaining arrays */
890a9d7e0ccSandi selinger   while (i<size1 && j<size2) {
891a9d7e0ccSandi selinger     if (in1[i] < in2[j]) {
892a9d7e0ccSandi selinger       out[l++] = in1[i++];
893a9d7e0ccSandi selinger     }
894a9d7e0ccSandi selinger     else if(in1[i] > in2[j]) {
895a9d7e0ccSandi selinger       out[l++] = in2[j++];
896a9d7e0ccSandi selinger     }
897a9d7e0ccSandi selinger     else {
898a9d7e0ccSandi selinger       out[l++] = in1[i];
899a9d7e0ccSandi selinger       i++, j++;
900a9d7e0ccSandi selinger     }
901a9d7e0ccSandi selinger   }
902a9d7e0ccSandi selinger 
903a9d7e0ccSandi selinger   while (i<size1 && k<size3) {
904a9d7e0ccSandi selinger     if (in1[i] < in3[k]) {
905a9d7e0ccSandi selinger       out[l++] = in1[i++];
906a9d7e0ccSandi selinger     }
907a9d7e0ccSandi selinger     else if(in1[i] > in3[k]) {
908a9d7e0ccSandi selinger       out[l++] = in3[k++];
909a9d7e0ccSandi selinger     }
910a9d7e0ccSandi selinger     else {
911a9d7e0ccSandi selinger       out[l++] = in1[i];
912a9d7e0ccSandi selinger       i++, k++;
913a9d7e0ccSandi selinger     }
914a9d7e0ccSandi selinger   }
915a9d7e0ccSandi selinger 
916a9d7e0ccSandi selinger   while (k<size3 && j<size2)  {
917a9d7e0ccSandi selinger     if (in3[k] < in2[j]) {
918a9d7e0ccSandi selinger       out[l++] = in3[k++];
919a9d7e0ccSandi selinger     }
920a9d7e0ccSandi selinger     else if(in3[k] > in2[j]) {
921a9d7e0ccSandi selinger       out[l++] = in2[j++];
922a9d7e0ccSandi selinger     }
923a9d7e0ccSandi selinger     else {
924a9d7e0ccSandi selinger       out[l++] = in3[k];
925a9d7e0ccSandi selinger       k++, j++;
926a9d7e0ccSandi selinger     }
927a9d7e0ccSandi selinger   }
928a9d7e0ccSandi selinger 
929a9d7e0ccSandi selinger   /* Traverse one remaining array */
930a9d7e0ccSandi selinger   while (i<size1) out[l++] = in1[i++];
931a9d7e0ccSandi selinger   while (j<size2) out[l++] = in2[j++];
932a9d7e0ccSandi selinger   while (k<size3) out[l++] = in3[k++];
933a9d7e0ccSandi selinger 
934a9d7e0ccSandi selinger   *size4 = l;
935a9d7e0ccSandi selinger }
936a9d7e0ccSandi selinger 
93773b67375Sandi selinger /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and  */
93873b67375Sandi selinger /* adds up the products. Two of these three multiplications are performed with existing (sequential)      */
93973b67375Sandi selinger /* matrix-matrix multiplications.  */
940a9d7e0ccSandi selinger #undef __FUNCT__
941a9d7e0ccSandi selinger #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI"
942a9d7e0ccSandi selinger PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C)
943a9d7e0ccSandi selinger {
944a9d7e0ccSandi selinger   PetscErrorCode     ierr;
945a9d7e0ccSandi selinger   MPI_Comm           comm;
946a9d7e0ccSandi selinger   PetscMPIInt        size;
947a9d7e0ccSandi selinger   Mat                Cmpi;
948a9d7e0ccSandi selinger   Mat_PtAPMPI        *ptap;
949a9d7e0ccSandi selinger   PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
950a9d7e0ccSandi selinger   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data;
951a9d7e0ccSandi selinger   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
952a9d7e0ccSandi selinger   Mat_MPIAIJ         *p        =(Mat_MPIAIJ*)P->data;
953a9d7e0ccSandi selinger   Mat_MPIAIJ         *c;
954a9d7e0ccSandi selinger   Mat_SeqAIJ         *adpd_seq, *p_off, *aopoth_seq;
955a9d7e0ccSandi selinger   PetscInt           adponz, adpdnz;
9567dbaa2c6Sandi selinger   PetscInt           *pi_loc,*dnz,*onz;
957a9d7e0ccSandi selinger   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
958a9d7e0ccSandi selinger   PetscInt           *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
959a9d7e0ccSandi selinger                      *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
960a016958eSandi selinger   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
961a9d7e0ccSandi selinger   PetscBT            lnkbt;
962a9d7e0ccSandi selinger   PetscScalar        *apa;
963a9d7e0ccSandi selinger   PetscReal          afill;
964a9d7e0ccSandi selinger   PetscMPIInt        rank;
965a9d7e0ccSandi selinger   Mat                adpd, aopoth;
966a9d7e0ccSandi selinger 
967a9d7e0ccSandi selinger   PetscFunctionBegin;
968a9d7e0ccSandi selinger   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
969a9d7e0ccSandi selinger   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
97073b67375Sandi selinger   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
971a016958eSandi selinger   ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend); CHKERRQ(ierr);
972a9d7e0ccSandi selinger 
973a9d7e0ccSandi selinger   /* create struct Mat_PtAPMPI and attached it to C later */
974a9d7e0ccSandi selinger   ierr = PetscNew(&ptap);CHKERRQ(ierr);
975a9d7e0ccSandi selinger 
976a9d7e0ccSandi selinger   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
977a9d7e0ccSandi selinger   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
978a9d7e0ccSandi selinger 
979a9d7e0ccSandi selinger   /* get P_loc by taking all local rows of P */
980a9d7e0ccSandi selinger   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
981a9d7e0ccSandi selinger 
982a9d7e0ccSandi selinger 
983a9d7e0ccSandi selinger   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
984a9d7e0ccSandi selinger   pi_loc = p_loc->i;
985a9d7e0ccSandi selinger 
986a9d7e0ccSandi selinger   /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
987a9d7e0ccSandi selinger   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
988a9d7e0ccSandi selinger   ierr      = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr);
989a9d7e0ccSandi selinger 
990a9d7e0ccSandi selinger   adpoi[0]    = 0;
991a9d7e0ccSandi selinger   ptap->api = api;
992a9d7e0ccSandi selinger   api[0] = 0;
993a9d7e0ccSandi selinger 
994a9d7e0ccSandi selinger   /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
995a9d7e0ccSandi selinger   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
996a9d7e0ccSandi selinger   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
997a9d7e0ccSandi selinger 
998a9d7e0ccSandi selinger   /* Symbolic calc of A_loc_diag * P_loc_diag */
99973b67375Sandi selinger   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);CHKERRQ(ierr);
1000a9d7e0ccSandi selinger   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1001a9d7e0ccSandi selinger   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1002a9d7e0ccSandi selinger   p_off = (Mat_SeqAIJ*)((p->B)->data);
1003a9d7e0ccSandi selinger   poff_i = p_off->i; poff_j = p_off->j;
1004a9d7e0ccSandi selinger 
100573b67375Sandi selinger   /* j_temp stores indices of a result row before they are added to the linked list */
1006a9d7e0ccSandi selinger   ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr);
1007a9d7e0ccSandi selinger 
1008a9d7e0ccSandi selinger 
1009a9d7e0ccSandi selinger   /* Symbolic calc of the A_diag * p_loc_off */
1010a9d7e0ccSandi selinger   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1011a9d7e0ccSandi selinger   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr);
1012a9d7e0ccSandi selinger   current_space = free_space_diag;
1013a9d7e0ccSandi selinger 
1014a9d7e0ccSandi selinger   for (i=0; i<am; i++) {
1015a9d7e0ccSandi selinger     /* A_diag * P_loc_off */
1016a9d7e0ccSandi selinger     nzi = adi[i+1] - adi[i];
1017a9d7e0ccSandi selinger     for (j=0; j<nzi; j++) {
1018a9d7e0ccSandi selinger       row  = *adj++;
1019a9d7e0ccSandi selinger       pnz  = poff_i[row+1] - poff_i[row];
1020a9d7e0ccSandi selinger       Jptr = poff_j + poff_i[row];
1021a9d7e0ccSandi selinger       for(i1 = 0; i1 < pnz; i1++) {
1022a9d7e0ccSandi selinger         j_temp[i1] = p->garray[Jptr[i1]];
1023a9d7e0ccSandi selinger       }
1024a9d7e0ccSandi selinger       /* add non-zero cols of P into the sorted linked list lnk */
1025a9d7e0ccSandi selinger       ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr);
1026a9d7e0ccSandi selinger     }
1027a9d7e0ccSandi selinger 
1028a9d7e0ccSandi selinger     adponz     = lnk[0];
1029a9d7e0ccSandi selinger     adpoi[i+1] = adpoi[i] + adponz;
1030a9d7e0ccSandi selinger 
1031a9d7e0ccSandi selinger     /* if free space is not available, double the total space in the list */
1032a9d7e0ccSandi selinger     if (current_space->local_remaining<adponz) {
1033a9d7e0ccSandi selinger       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1034a9d7e0ccSandi selinger       nspacedouble++;
1035a9d7e0ccSandi selinger     }
1036a9d7e0ccSandi selinger 
1037a9d7e0ccSandi selinger     /* Copy data into free space, then initialize lnk */
1038a9d7e0ccSandi selinger     ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1039a9d7e0ccSandi selinger 
1040a9d7e0ccSandi selinger     current_space->array           += adponz;
1041a9d7e0ccSandi selinger     current_space->local_used      += adponz;
1042a9d7e0ccSandi selinger     current_space->local_remaining -= adponz;
1043a9d7e0ccSandi selinger   }
1044a9d7e0ccSandi selinger 
1045a9d7e0ccSandi selinger   /* Symbolic calc of A_off * P_oth */
104673b67375Sandi selinger   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);CHKERRQ(ierr);
1047a9d7e0ccSandi selinger   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1048a9d7e0ccSandi selinger   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;
1049a9d7e0ccSandi selinger 
1050a9d7e0ccSandi selinger   /* Allocate space for apj, adpj, aopj, ... */
1051a9d7e0ccSandi selinger   /* destroy lists of free space and other temporary array(s) */
1052a9d7e0ccSandi selinger 
1053a9d7e0ccSandi selinger   ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr);
1054a9d7e0ccSandi selinger   ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr);
1055a9d7e0ccSandi selinger 
1056a9d7e0ccSandi selinger   /* Copy from linked list to j-array */
1057a9d7e0ccSandi selinger   ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr);
1058a9d7e0ccSandi selinger   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1059a9d7e0ccSandi selinger 
1060a9d7e0ccSandi selinger   adpoJ = adpoj;
1061a9d7e0ccSandi selinger   adpdJ = adpdj;
1062a9d7e0ccSandi selinger   aopJ = aopothj;
1063a9d7e0ccSandi selinger   apj  = ptap->apj;
1064a9d7e0ccSandi selinger   apJ = apj; /* still empty */
1065a9d7e0ccSandi selinger 
1066a9d7e0ccSandi selinger   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1067a9d7e0ccSandi selinger   /* A_diag * P_loc_diag to get A*P */
1068a9d7e0ccSandi selinger   for (i = 0; i < am; i++) {
1069a9d7e0ccSandi selinger     aopnz  =  aopothi[i+1] -  aopothi[i];
1070a9d7e0ccSandi selinger     adponz = adpoi[i+1] - adpoi[i];
1071a9d7e0ccSandi selinger     adpdnz = adpdi[i+1] - adpdi[i];
1072a9d7e0ccSandi selinger 
1073a9d7e0ccSandi selinger     /* Correct indices from A_diag*P_diag */
1074a9d7e0ccSandi selinger     for(i1 = 0; i1 < adpdnz; i1++) {
1075a016958eSandi selinger       adpdJ[i1] += p_colstart;
1076a9d7e0ccSandi selinger     }
1077a9d7e0ccSandi selinger     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1078a9d7e0ccSandi selinger     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1079a9d7e0ccSandi selinger     ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz); CHKERRQ(ierr);
1080a9d7e0ccSandi selinger 
1081a9d7e0ccSandi selinger     aopJ += aopnz;
1082a9d7e0ccSandi selinger     adpoJ += adponz;
1083a9d7e0ccSandi selinger     adpdJ += adpdnz;
1084a9d7e0ccSandi selinger     apJ += apnz;
1085a9d7e0ccSandi selinger     api[i+1] = api[i] + apnz;
1086a9d7e0ccSandi selinger   }
1087a9d7e0ccSandi selinger 
1088a9d7e0ccSandi selinger   /* malloc apa to store dense row A[i,:]*P */
1089a9d7e0ccSandi selinger   ierr = PetscCalloc1(pN+2,&apa);CHKERRQ(ierr);
1090a9d7e0ccSandi selinger 
1091a9d7e0ccSandi selinger   ptap->apa = apa;
1092a9d7e0ccSandi selinger   /* create and assemble symbolic parallel matrix Cmpi */
1093a9d7e0ccSandi selinger   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1094a9d7e0ccSandi selinger   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1095a9d7e0ccSandi selinger   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
1096a9d7e0ccSandi selinger 
1097a9d7e0ccSandi selinger   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1098a9d7e0ccSandi selinger   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
109917f75172Sandi selinger 
110017f75172Sandi selinger 
1101*904d1e70Sandi selinger   ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr);
1102a9d7e0ccSandi selinger   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1103a9d7e0ccSandi selinger   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11047dbaa2c6Sandi selinger   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1105a9d7e0ccSandi selinger 
110617f75172Sandi selinger 
1107a9d7e0ccSandi selinger   ptap->destroy        = Cmpi->ops->destroy;
1108a9d7e0ccSandi selinger   ptap->duplicate      = Cmpi->ops->duplicate;
1109a9d7e0ccSandi selinger   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1110a9d7e0ccSandi selinger   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
1111a9d7e0ccSandi selinger   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
1112a9d7e0ccSandi selinger 
1113a9d7e0ccSandi selinger   /* attach the supporting struct to Cmpi for reuse */
1114a9d7e0ccSandi selinger   c       = (Mat_MPIAIJ*)Cmpi->data;
1115a9d7e0ccSandi selinger   c->ptap = ptap;
1116a9d7e0ccSandi selinger   *C = Cmpi;
1117a9d7e0ccSandi selinger 
1118a9d7e0ccSandi selinger   /* set MatInfo */
111973b67375Sandi selinger   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
112073b67375Sandi selinger   if (afill < 1.0) afill = 1.0;
1121a9d7e0ccSandi selinger   Cmpi->info.mallocs           = nspacedouble;
1122a9d7e0ccSandi selinger   Cmpi->info.fill_ratio_given  = fill;
1123a9d7e0ccSandi selinger   Cmpi->info.fill_ratio_needed = afill;
1124a9d7e0ccSandi selinger 
1125a9d7e0ccSandi selinger #if defined(PETSC_USE_INFO)
1126a9d7e0ccSandi selinger   if (api[am]) {
1127a9d7e0ccSandi selinger     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1128a9d7e0ccSandi selinger     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1129a9d7e0ccSandi selinger   } else {
1130a9d7e0ccSandi selinger     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1131a9d7e0ccSandi selinger   }
1132a9d7e0ccSandi selinger #endif
1133a9d7e0ccSandi selinger 
113473b67375Sandi selinger   ierr = MatDestroy(&aopoth);CHKERRQ(ierr);
113573b67375Sandi selinger   ierr = MatDestroy(&adpd);CHKERRQ(ierr);
113673b67375Sandi selinger   ierr = PetscFree(j_temp);CHKERRQ(ierr);
113773b67375Sandi selinger   ierr = PetscFree(adpoj);CHKERRQ(ierr);
113873b67375Sandi selinger   ierr = PetscFree(adpoi);CHKERRQ(ierr);
1139a9d7e0ccSandi selinger   PetscFunctionReturn(0);
1140a9d7e0ccSandi selinger }
1141a9d7e0ccSandi selinger 
1142a9d7e0ccSandi selinger 
1143f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
1144c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1145f96d37fcSHong Zhang {
1146f96d37fcSHong Zhang   PetscErrorCode ierr;
1147c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
11484ede91e3SHong Zhang   PetscInt       aN=A->cmap->N,alg=1; /* set default algorithm */
11494ede91e3SHong Zhang   PetscBool      flg;
1150f96d37fcSHong Zhang 
1151f96d37fcSHong Zhang   PetscFunctionBegin;
1152c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
1153c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
115468eacb73SHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
11554ede91e3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);CHKERRQ(ierr);
1156c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1157c216dbf3SHong Zhang 
1158669e9c86SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
1159c216dbf3SHong Zhang     switch (alg) {
1160c216dbf3SHong Zhang     case 1:
11614ede91e3SHong Zhang       if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */
11624ede91e3SHong Zhang         MatInfo     Ainfo,Pinfo;
11634ede91e3SHong Zhang         PetscInt    nz_local;
11644ede91e3SHong Zhang         PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
11654ede91e3SHong Zhang         MPI_Comm    comm;
11664ede91e3SHong Zhang 
11674ede91e3SHong Zhang         ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
11684ede91e3SHong Zhang         ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
11694ede91e3SHong Zhang         nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */
11704ede91e3SHong Zhang 
11714ede91e3SHong Zhang         if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
11724ede91e3SHong Zhang         ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
11734ede91e3SHong Zhang         ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
11744ede91e3SHong Zhang 
11754ede91e3SHong Zhang         if (alg_scalable) {
11764ede91e3SHong Zhang           alg  = 0; /* scalable algorithm would slower than nonscalable algorithm */
11774ede91e3SHong Zhang           ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
11784ede91e3SHong Zhang           break;
11794ede91e3SHong Zhang         }
11804ede91e3SHong Zhang       }
11812bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
1182c216dbf3SHong Zhang       break;
1183c216dbf3SHong Zhang     case 2:
1184275476c6SMatthew G. Knepley     {
1185ab1b013aSHong Zhang       Mat         Pt;
1186ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
1187ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
1188ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
1189ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
1190ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
1191ab1b013aSHong Zhang       ptap     = c->ptap;
1192ab1b013aSHong Zhang       ptap->Pt = Pt;
1193c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
1194c216dbf3SHong Zhang       PetscFunctionReturn(0);
1195275476c6SMatthew G. Knepley     }
1196c216dbf3SHong Zhang       break;
11974ede91e3SHong Zhang     default: /* scalable algorithm */
11986da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
1199c216dbf3SHong Zhang       break;
1200c216dbf3SHong Zhang     }
1201669e9c86SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
1202c216dbf3SHong Zhang   }
1203669e9c86SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
1204c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
1205669e9c86SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
1206ab1b013aSHong Zhang   PetscFunctionReturn(0);
1207ab1b013aSHong Zhang }
1208ab1b013aSHong Zhang 
1209c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
1210c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1211c216dbf3SHong Zhang {
1212c216dbf3SHong Zhang   PetscErrorCode ierr;
12132bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
12142bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
12152bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
1216c216dbf3SHong Zhang 
1217c216dbf3SHong Zhang   PetscFunctionBegin;
1218c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
1219c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
1220f96d37fcSHong Zhang   PetscFunctionReturn(0);
1221f96d37fcSHong Zhang }
1222f96d37fcSHong Zhang 
1223c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat,MatDuplicateOption,Mat*);
12244e201bd7SHong Zhang 
1225c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1226912837bbSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1227912837bbSHong Zhang {
1228912837bbSHong Zhang   PetscErrorCode      ierr;
1229912837bbSHong Zhang   Mat_PtAPMPI         *ptap;
12307f3e1f30SHong Zhang   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1231912837bbSHong Zhang   MPI_Comm            comm;
1232912837bbSHong Zhang   PetscMPIInt         size,rank;
1233912837bbSHong Zhang   Mat                 Cmpi;
1234912837bbSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
12354ede91e3SHong Zhang   PetscInt            pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
12367f3e1f30SHong Zhang   PetscInt            *lnk,i,k,nsend;
1237912837bbSHong Zhang   PetscBT             lnkbt;
1238912837bbSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1239912837bbSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
12407f3e1f30SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi;
1241912837bbSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1242912837bbSHong Zhang   MPI_Request         *swaits,*rwaits;
1243912837bbSHong Zhang   MPI_Status          *sstatus,rstatus;
1244912837bbSHong Zhang   PetscLayout         rowmap;
1245912837bbSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1246912837bbSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
12477f3e1f30SHong Zhang   PetscInt            *Jptr,*prmap=p->garray,con,j,Crmax;
12484ede91e3SHong Zhang   Mat_SeqAIJ          *a_loc,*c_loc,*c_oth;
1249912837bbSHong Zhang   PetscTable          ta;
1250912837bbSHong Zhang 
1251912837bbSHong Zhang   PetscFunctionBegin;
1252912837bbSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1253912837bbSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1254912837bbSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1255912837bbSHong Zhang 
1256912837bbSHong Zhang   /* create symbolic parallel matrix Cmpi */
1257912837bbSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1258912837bbSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1259912837bbSHong Zhang 
1260912837bbSHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1261912837bbSHong Zhang 
1262912837bbSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1263912837bbSHong Zhang   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1264912837bbSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
1265912837bbSHong Zhang 
1266912837bbSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1267912837bbSHong Zhang   /* --------------------------------- */
1268912837bbSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1269912837bbSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1270912837bbSHong Zhang 
1271669e9c86SHong Zhang   /* (1) compute symbolic A_loc */
1272669e9c86SHong Zhang   /* ---------------------------*/
12734ede91e3SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr);
1274912837bbSHong Zhang 
12754ede91e3SHong Zhang   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1276912837bbSHong Zhang   /* ------------------------------------ */
12774ede91e3SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
1278912837bbSHong Zhang 
1279912837bbSHong Zhang   /* (3) send coj of C_oth to other processors  */
1280912837bbSHong Zhang   /* ------------------------------------------ */
1281912837bbSHong Zhang   /* determine row ownership */
1282912837bbSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
1283912837bbSHong Zhang   rowmap->n  = pn;
1284912837bbSHong Zhang   rowmap->bs = 1;
1285912837bbSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
1286912837bbSHong Zhang   owners = rowmap->range;
1287912837bbSHong Zhang 
1288912837bbSHong Zhang   /* determine the number of messages to send, their lengths */
1289912837bbSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
1290912837bbSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1291912837bbSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1292912837bbSHong Zhang 
1293912837bbSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1294912837bbSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
1295912837bbSHong Zhang   con   = ptap->C_oth->rmap->n;
1296912837bbSHong Zhang   proc  = 0;
1297912837bbSHong Zhang   for (i=0; i<con; i++) {
1298912837bbSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
12997f3e1f30SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*A) to be sent to [proc] */
1300912837bbSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1301912837bbSHong Zhang   }
1302912837bbSHong Zhang 
1303912837bbSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
1304912837bbSHong Zhang   owners_co[0] = 0;
1305912837bbSHong Zhang   nsend        = 0;
1306912837bbSHong Zhang   for (proc=0; proc<size; proc++) {
1307912837bbSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1308912837bbSHong Zhang     if (len_s[proc]) {
1309912837bbSHong Zhang       nsend++;
1310912837bbSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1311912837bbSHong Zhang       len         += len_si[proc];
1312912837bbSHong Zhang     }
1313912837bbSHong Zhang   }
1314912837bbSHong Zhang 
1315912837bbSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
1316912837bbSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
1317912837bbSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
1318912837bbSHong Zhang 
1319912837bbSHong Zhang   /* post the Irecv and Isend of coj */
1320912837bbSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1321912837bbSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1322912837bbSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
1323912837bbSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1324912837bbSHong Zhang     if (!len_s[proc]) continue;
1325912837bbSHong Zhang     i    = owners_co[proc];
1326912837bbSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1327912837bbSHong Zhang     k++;
1328912837bbSHong Zhang   }
1329912837bbSHong Zhang 
13307f3e1f30SHong Zhang   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1331912837bbSHong Zhang   /* ---------------------------------------- */
13324ede91e3SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
1333912837bbSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1334912837bbSHong Zhang 
1335912837bbSHong Zhang   /* receives coj are complete */
1336912837bbSHong Zhang   for (i=0; i<nrecv; i++) {
1337912837bbSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1338912837bbSHong Zhang   }
1339912837bbSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1340912837bbSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1341912837bbSHong Zhang 
1342912837bbSHong Zhang   /* add received column indices into ta to update Crmax */
13434ede91e3SHong Zhang   a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data;
1344669e9c86SHong Zhang 
1345669e9c86SHong Zhang   /* create and initialize a linked list */
13464ede91e3SHong Zhang   ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */
13474ede91e3SHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta);
1348669e9c86SHong Zhang 
1349912837bbSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
1350912837bbSHong Zhang     Jptr = buf_rj[k];
1351912837bbSHong Zhang     for (j=0; j<len_r[k]; j++) {
1352912837bbSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
1353912837bbSHong Zhang     }
1354912837bbSHong Zhang   }
1355912837bbSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
1356912837bbSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
1357912837bbSHong Zhang 
1358912837bbSHong Zhang   /* (4) send and recv coi */
1359912837bbSHong Zhang   /*-----------------------*/
1360912837bbSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1361912837bbSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1362912837bbSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1363912837bbSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1364912837bbSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1365912837bbSHong Zhang     if (!len_s[proc]) continue;
1366912837bbSHong Zhang     /* form outgoing message for i-structure:
1367912837bbSHong Zhang          buf_si[0]:                 nrows to be sent
1368912837bbSHong Zhang                [1:nrows]:           row index (global)
1369912837bbSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1370912837bbSHong Zhang     */
1371912837bbSHong Zhang     /*-------------------------------------------*/
1372912837bbSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1373912837bbSHong Zhang     buf_si_i    = buf_si + nrows+1;
1374912837bbSHong Zhang     buf_si[0]   = nrows;
1375912837bbSHong Zhang     buf_si_i[0] = 0;
1376912837bbSHong Zhang     nrows       = 0;
1377912837bbSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1378912837bbSHong Zhang       nzi = coi[i+1] - coi[i];
1379912837bbSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1380912837bbSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1381912837bbSHong Zhang       nrows++;
1382912837bbSHong Zhang     }
1383912837bbSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1384912837bbSHong Zhang     k++;
1385912837bbSHong Zhang     buf_si += len_si[proc];
1386912837bbSHong Zhang   }
1387912837bbSHong Zhang   for (i=0; i<nrecv; i++) {
1388912837bbSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1389912837bbSHong Zhang   }
1390912837bbSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1391912837bbSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
1392912837bbSHong Zhang 
1393912837bbSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
1394912837bbSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1395912837bbSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1396912837bbSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1397912837bbSHong Zhang 
1398912837bbSHong Zhang   /* (5) compute the local portion of Cmpi      */
1399912837bbSHong Zhang   /* ------------------------------------------ */
1400912837bbSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1401912837bbSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
1402912837bbSHong Zhang   current_space = free_space;
1403912837bbSHong Zhang 
1404912837bbSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
1405912837bbSHong Zhang   for (k=0; k<nrecv; k++) {
1406912837bbSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1407912837bbSHong Zhang     nrows       = *buf_ri_k[k];
1408912837bbSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1409912837bbSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1410912837bbSHong Zhang   }
1411912837bbSHong Zhang 
14124ede91e3SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr);
14134ede91e3SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1414912837bbSHong Zhang   for (i=0; i<pn; i++) {
1415912837bbSHong Zhang     /* add C_loc into Cmpi */
1416912837bbSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
1417912837bbSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
1418912837bbSHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1419912837bbSHong Zhang 
1420912837bbSHong Zhang     /* add received col data into lnk */
1421912837bbSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
1422912837bbSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1423912837bbSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1424912837bbSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1425912837bbSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1426912837bbSHong Zhang         nextrow[k]++; nextci[k]++;
1427912837bbSHong Zhang       }
1428912837bbSHong Zhang     }
1429912837bbSHong Zhang     nzi = lnk[0];
1430912837bbSHong Zhang 
1431912837bbSHong Zhang     /* copy data into free space, then initialize lnk */
14324ede91e3SHong Zhang     ierr = PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1433912837bbSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
1434912837bbSHong Zhang   }
1435912837bbSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1436912837bbSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1437912837bbSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
1438912837bbSHong Zhang 
1439912837bbSHong Zhang   /* local sizes and preallocation */
14404ede91e3SHong Zhang   ierr = MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1441912837bbSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1442912837bbSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1443912837bbSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1444912837bbSHong Zhang 
1445912837bbSHong Zhang   /* members in merge */
1446912837bbSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
1447912837bbSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
1448912837bbSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1449912837bbSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1450912837bbSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1451912837bbSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1452912837bbSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
1453912837bbSHong Zhang 
1454912837bbSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1455912837bbSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
1456912837bbSHong Zhang   c->ptap         = ptap;
1457912837bbSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1458912837bbSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1459912837bbSHong Zhang 
1460912837bbSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1461912837bbSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
1462912837bbSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1463669e9c86SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
1464912837bbSHong Zhang   *C                     = Cmpi;
1465912837bbSHong Zhang   PetscFunctionReturn(0);
1466912837bbSHong Zhang }
1467912837bbSHong Zhang 
1468912837bbSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1469912837bbSHong Zhang {
1470912837bbSHong Zhang   PetscErrorCode    ierr;
1471669e9c86SHong Zhang   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1472669e9c86SHong Zhang   Mat_SeqAIJ        *c_seq;
1473912837bbSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
1474669e9c86SHong Zhang   Mat               A_loc,C_loc,C_oth;
1475912837bbSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
1476912837bbSHong Zhang   const PetscInt    *cols;
1477912837bbSHong Zhang   const PetscScalar *vals;
1478912837bbSHong Zhang 
1479912837bbSHong Zhang   PetscFunctionBegin;
1480912837bbSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1481912837bbSHong Zhang 
1482912837bbSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
1483669e9c86SHong Zhang     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1484669e9c86SHong Zhang     /* 1) get R = Pd^T, Ro = Po^T */
1485669e9c86SHong Zhang     /*----------------------------*/
1486912837bbSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1487912837bbSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1488669e9c86SHong Zhang 
1489669e9c86SHong Zhang     /* 2) compute numeric A_loc */
1490669e9c86SHong Zhang     /*--------------------------*/
14914ede91e3SHong Zhang     ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr);
1492912837bbSHong Zhang   }
1493912837bbSHong Zhang 
1494669e9c86SHong Zhang   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
14954ede91e3SHong Zhang   A_loc = ptap->A_loc;
1496669e9c86SHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr);
1497669e9c86SHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr);
1498912837bbSHong Zhang   C_loc = ptap->C_loc;
1499912837bbSHong Zhang   C_oth = ptap->C_oth;
1500912837bbSHong Zhang 
1501912837bbSHong Zhang   /* add C_loc and Co to to C */
1502912837bbSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
1503912837bbSHong Zhang 
1504912837bbSHong Zhang   /* C_loc -> C */
1505912837bbSHong Zhang   cm    = C_loc->rmap->N;
1506912837bbSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
1507912837bbSHong Zhang   cols = c_seq->j;
1508912837bbSHong Zhang   vals = c_seq->a;
1509912837bbSHong Zhang   for (i=0; i<cm; i++) {
1510912837bbSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
1511912837bbSHong Zhang     row = rstart + i;
1512912837bbSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1513912837bbSHong Zhang     cols += ncols; vals += ncols;
1514912837bbSHong Zhang   }
1515912837bbSHong Zhang 
1516912837bbSHong Zhang   /* Co -> C, off-processor part */
1517912837bbSHong Zhang   cm    = C_oth->rmap->N;
1518912837bbSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
1519912837bbSHong Zhang   cols  = c_seq->j;
1520912837bbSHong Zhang   vals  = c_seq->a;
1521912837bbSHong Zhang   for (i=0; i<cm; i++) {
1522912837bbSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
1523912837bbSHong Zhang     row = p->garray[i];
1524912837bbSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
1525912837bbSHong Zhang     cols += ncols; vals += ncols;
1526912837bbSHong Zhang   }
1527912837bbSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1528912837bbSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1529912837bbSHong Zhang 
1530912837bbSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
1531912837bbSHong Zhang   PetscFunctionReturn(0);
1532912837bbSHong Zhang }
1533912837bbSHong Zhang 
15346da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1535d6ab1dc0SHong Zhang {
1536d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1537d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1538d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1539d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1540d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1541e745005bSHong Zhang   PetscInt            *adj;
1542e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1543e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1544d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1545ce94432eSBarry Smith   MPI_Comm            comm;
1546d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1547d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1548d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1549d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1550d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1551d6ab1dc0SHong Zhang   MPI_Status          *status;
1552d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
15534e201bd7SHong Zhang   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1554d6ab1dc0SHong Zhang   Mat                 A_loc;
1555d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1556d6ab1dc0SHong Zhang 
1557d6ab1dc0SHong Zhang   PetscFunctionBegin;
1558ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1559d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1560d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1561d6ab1dc0SHong Zhang 
1562d6ab1dc0SHong Zhang   ptap  = c->ptap;
1563d6ab1dc0SHong Zhang   merge = ptap->merge;
1564d6ab1dc0SHong Zhang 
1565e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1566e745005bSHong Zhang   /*------------------------------------------*/
1567d6ab1dc0SHong Zhang   /* get data from symbolic products */
1568d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1569854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1570d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1571d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
15721795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1573d6ab1dc0SHong Zhang 
1574d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1575d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1576d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1577d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1578d6ab1dc0SHong Zhang   ai    = a_loc->i;
1579d6ab1dc0SHong Zhang   aj    = a_loc->j;
1580d6ab1dc0SHong Zhang 
1581d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1582d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1583d6ab1dc0SHong Zhang     adj = aj + ai[i];
1584d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1585d6ab1dc0SHong Zhang 
1586d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1587e745005bSHong Zhang     /*-------------------------------------------------------------*/
1588d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1589d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1590d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1591d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1592d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1593d6ab1dc0SHong Zhang       row = poJ[j];
1594d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1595d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1596e745005bSHong Zhang       /* perform sparse axpy */
1597e745005bSHong Zhang       nexta  = 0;
1598d6ab1dc0SHong Zhang       valtmp = pA[j];
1599e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1600e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1601e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1602e745005bSHong Zhang           nexta++;
1603d6ab1dc0SHong Zhang         }
1604e745005bSHong Zhang       }
1605e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1606d6ab1dc0SHong Zhang     }
1607d6ab1dc0SHong Zhang 
1608d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1609d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1610d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1611d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1612d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1613d6ab1dc0SHong Zhang       row = pdJ[j];
1614d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1615d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1616e745005bSHong Zhang       /* perform sparse axpy */
1617e745005bSHong Zhang       nexta  = 0;
1618d6ab1dc0SHong Zhang       valtmp = pA[j];
1619e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1620e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1621e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1622e745005bSHong Zhang           nexta++;
1623d6ab1dc0SHong Zhang         }
1624e745005bSHong Zhang       }
1625e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1626d6ab1dc0SHong Zhang     }
1627d6ab1dc0SHong Zhang   }
1628d6ab1dc0SHong Zhang 
1629d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1630d6ab1dc0SHong Zhang   /*------------------------------------*/
1631d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1632d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1633d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1634d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1635d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1636d6ab1dc0SHong Zhang 
1637dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1638d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1639d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1640d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1641d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1642d6ab1dc0SHong Zhang     k++;
1643d6ab1dc0SHong Zhang   }
1644d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1645d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1646d6ab1dc0SHong Zhang 
1647d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1648d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1649d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1650d6ab1dc0SHong Zhang 
1651d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1652d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1653dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1654d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1655e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1656d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1657d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1658d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1659d6ab1dc0SHong Zhang   }
1660d6ab1dc0SHong Zhang 
1661d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1662d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1663d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1664d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1665d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1666d6ab1dc0SHong Zhang     /* add received vals into ba */
1667d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1668d6ab1dc0SHong Zhang       /* i-th row */
1669d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1670d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1671d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1672d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1673d6ab1dc0SHong Zhang         nextcj = 0;
1674d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1675d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1676d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1677d6ab1dc0SHong Zhang           }
1678d6ab1dc0SHong Zhang         }
1679d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1680d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1681d6ab1dc0SHong Zhang       }
1682d6ab1dc0SHong Zhang     }
1683d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1684d6ab1dc0SHong Zhang   }
1685d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1686d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1687d6ab1dc0SHong Zhang 
1688d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1689d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1690d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1691d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1692d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1693d6ab1dc0SHong Zhang }
1694d6ab1dc0SHong Zhang 
16956da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1696d6ab1dc0SHong Zhang {
1697d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1698d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1699d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
17000298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
170133ba5e2fSHong Zhang   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1702d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1703d6ab1dc0SHong Zhang   PetscInt            nnz;
1704d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1705d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1706ce94432eSBarry Smith   MPI_Comm            comm;
1707d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1708d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1709d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1710d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1711d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1712d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1713d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1714d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1715d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1716d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
17174b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1718d6ab1dc0SHong Zhang   PetscScalar         *vals;
1719d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
17200acc65b4SHong Zhang   PetscTable          ta;
1721d6ab1dc0SHong Zhang 
1722d6ab1dc0SHong Zhang   PetscFunctionBegin;
1723ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1724d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
17256c4ed002SBarry 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);
1726d6ab1dc0SHong Zhang 
1727d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1728d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1729d6ab1dc0SHong Zhang 
1730d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1731b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1732d6ab1dc0SHong Zhang 
1733d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1734d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
17352205254eSKarl Rupp 
1736d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1737d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1738d6ab1dc0SHong Zhang   ai          = a_loc->i;
1739d6ab1dc0SHong Zhang   aj          = a_loc->j;
1740d6ab1dc0SHong Zhang 
1741d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1742d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1743d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1745d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1746d6ab1dc0SHong Zhang 
1747d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1748d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1749d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1750d6ab1dc0SHong Zhang 
1751d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1752d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1753d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1754854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1755d6ab1dc0SHong Zhang   coi[0] = 0;
1756d6ab1dc0SHong Zhang 
1757d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1758f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1759a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1760d6ab1dc0SHong Zhang   current_space = free_space;
176119f0e57aSHong Zhang 
1762d6ab1dc0SHong Zhang   /* create and initialize a linked list */
176333ba5e2fSHong Zhang   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
17644b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
17650acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
176633ba5e2fSHong Zhang 
17670acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1768d6ab1dc0SHong Zhang 
1769d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1770d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1771d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1772d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1773d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1774d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1775d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1776d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1777d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1778d6ab1dc0SHong Zhang     }
1779d6ab1dc0SHong Zhang     nnz = lnk[0];
1780d6ab1dc0SHong Zhang 
1781d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1782d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1783f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1784d6ab1dc0SHong Zhang       nspacedouble++;
1785d6ab1dc0SHong Zhang     }
1786d6ab1dc0SHong Zhang 
1787d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1788d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
17892205254eSKarl Rupp 
1790d6ab1dc0SHong Zhang     current_space->array           += nnz;
1791d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1792d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
17932205254eSKarl Rupp 
1794d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1795d6ab1dc0SHong Zhang   }
1796d6ab1dc0SHong Zhang 
1797854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1798d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17990acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
18002205254eSKarl Rupp 
1801118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1802d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1803d6ab1dc0SHong Zhang 
1804d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1805d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1806d6ab1dc0SHong Zhang   /* determine row ownership */
1807b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1808d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
18092205254eSKarl Rupp 
1810d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1811d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
18122205254eSKarl Rupp 
1813d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1815d6ab1dc0SHong Zhang 
1816d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
18171795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1818785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
18192205254eSKarl Rupp 
1820d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1821d6ab1dc0SHong Zhang   merge->nsend = 0;
1822d6ab1dc0SHong Zhang 
1823854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1824d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1825d6ab1dc0SHong Zhang 
1826d6ab1dc0SHong Zhang   proc = 0;
1827d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1828d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1829d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1830d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1831d6ab1dc0SHong Zhang   }
1832d6ab1dc0SHong Zhang 
1833d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1834d6ab1dc0SHong Zhang   owners_co[0] = 0;
1835d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1836d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1837d6ab1dc0SHong Zhang     if (len_si[proc]) {
1838d6ab1dc0SHong Zhang       merge->nsend++;
1839d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1840d6ab1dc0SHong Zhang       len         += len_si[proc];
1841d6ab1dc0SHong Zhang     }
1842d6ab1dc0SHong Zhang   }
1843d6ab1dc0SHong Zhang 
1844d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
18450298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1846d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1847d6ab1dc0SHong Zhang 
1848d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1849d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1850d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1851854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1852d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1853d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1854d6ab1dc0SHong Zhang     i    = owners_co[proc];
1855d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1856d6ab1dc0SHong Zhang     k++;
1857d6ab1dc0SHong Zhang   }
1858d6ab1dc0SHong Zhang 
1859d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1860785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1861d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1862c280ed6aSJed Brown     PetscMPIInt icompleted;
1863c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1864d6ab1dc0SHong Zhang   }
1865d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1866d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1867d6ab1dc0SHong Zhang 
18680acc65b4SHong Zhang   /* add received column indices into table to update Armax */
186933ba5e2fSHong Zhang   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
18700acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
18710acc65b4SHong Zhang     Jptr = buf_rj[k];
18720acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
18730acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
18740acc65b4SHong Zhang     }
18750acc65b4SHong Zhang   }
18760acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
187733ba5e2fSHong Zhang   /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */
18780acc65b4SHong Zhang 
1879d6ab1dc0SHong Zhang   /* send and recv coi */
1880d6ab1dc0SHong Zhang   /*-------------------*/
1881d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1882d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1883854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1884d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1885d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1886d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1887d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1888d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1889d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1890d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1891d6ab1dc0SHong Zhang     */
1892d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1893d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1894d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1895d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1896d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1897d6ab1dc0SHong Zhang     nrows       = 0;
1898d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1899d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1900d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1901d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1902d6ab1dc0SHong Zhang       nrows++;
1903d6ab1dc0SHong Zhang     }
1904d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1905d6ab1dc0SHong Zhang     k++;
1906d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1907d6ab1dc0SHong Zhang   }
1908d6ab1dc0SHong Zhang   i = merge->nrecv;
1909d6ab1dc0SHong Zhang   while (i--) {
1910c280ed6aSJed Brown     PetscMPIInt icompleted;
1911c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1912d6ab1dc0SHong Zhang   }
1913d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1914d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1915d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1916d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1917d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1918d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1919d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1920d6ab1dc0SHong Zhang 
1921d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1922d6ab1dc0SHong Zhang   /*------------------------------------------*/
1923d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1924854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1925d6ab1dc0SHong Zhang   bi[0] = 0;
1926d6ab1dc0SHong Zhang 
1927d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1928f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1929a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1930d6ab1dc0SHong Zhang   current_space = free_space;
1931d6ab1dc0SHong Zhang 
1932dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1933d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1934d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1935d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1936d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
19372205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1938d6ab1dc0SHong Zhang   }
1939d6ab1dc0SHong Zhang 
19400acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1941d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1942d6ab1dc0SHong Zhang   rmax = 0;
1943d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1944d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1945d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1946d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1947d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1948d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1949d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1950d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1951d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1952d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1953d6ab1dc0SHong Zhang     }
1954d6ab1dc0SHong Zhang 
1955d6ab1dc0SHong Zhang     /* add received col data into lnk */
1956d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1957d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1958d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1959d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1960d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1961d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1962d6ab1dc0SHong Zhang       }
1963d6ab1dc0SHong Zhang     }
1964d6ab1dc0SHong Zhang     nnz = lnk[0];
1965d6ab1dc0SHong Zhang 
1966d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1967d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1968f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1969d6ab1dc0SHong Zhang       nspacedouble++;
1970d6ab1dc0SHong Zhang     }
1971d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1972d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1973d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
19742205254eSKarl Rupp 
1975d6ab1dc0SHong Zhang     current_space->array           += nnz;
1976d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1977d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
19782205254eSKarl Rupp 
1979d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1980d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1981d6ab1dc0SHong Zhang   }
1982f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1983d6ab1dc0SHong Zhang 
1984854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1985d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1986118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1987d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1988d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
19890acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
19900acc65b4SHong Zhang 
1991d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1992d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1993d6ab1dc0SHong Zhang 
1994d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1995d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
19961795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1997d6ab1dc0SHong Zhang 
1998d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1999d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
200033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
2001d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
2002d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
2003d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2004d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
2005d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
2006d6ab1dc0SHong Zhang     row  = i + rstart;
2007d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
2008d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
2009d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
2010d6ab1dc0SHong Zhang   }
2011d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2012d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2013d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
2014d6ab1dc0SHong Zhang 
2015d6ab1dc0SHong Zhang   merge->bi        = bi;
2016d6ab1dc0SHong Zhang   merge->bj        = bj;
2017d6ab1dc0SHong Zhang   merge->coi       = coi;
2018d6ab1dc0SHong Zhang   merge->coj       = coj;
2019d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
2020d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
2021d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
2022d6ab1dc0SHong Zhang 
2023d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
2024d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
20252205254eSKarl Rupp 
2026d6ab1dc0SHong Zhang   c->ptap     = ptap;
20270298fd71SBarry Smith   ptap->api   = NULL;
20280298fd71SBarry Smith   ptap->apj   = NULL;
2029d6ab1dc0SHong Zhang   ptap->merge = merge;
20300298fd71SBarry Smith   ptap->apa   = NULL;
2031a560ef98SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
2032a560ef98SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
2033a560ef98SHong Zhang 
2034a560ef98SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2035a560ef98SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
2036a560ef98SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
2037f96d37fcSHong Zhang 
2038f96d37fcSHong Zhang   *C = Cmpi;
2039f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
2040f96d37fcSHong Zhang   if (bi[pn] != 0) {
204157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
204257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
2043f96d37fcSHong Zhang   } else {
2044f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
2045f96d37fcSHong Zhang   }
2046f96d37fcSHong Zhang #endif
2047f96d37fcSHong Zhang   PetscFunctionReturn(0);
2048f96d37fcSHong Zhang }
2049