1*06ef35abSSatish Balay /*$Id: baijov.c,v 1.61 2001/03/09 19:03:09 balay Exp balay $*/ 2ca161407SBarry Smith 3d5b485abSSatish Balay /* 4d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 5d5b485abSSatish Balay and to find submatrices that were shared across processors. 6d5b485abSSatish Balay */ 770f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h" 80a835dfdSSatish Balay #include "petscbt.h" 9d5b485abSSatish Balay 10df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat,int,IS *); 11df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat,int,char **,int*,int**); 12df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat,int,int **,int**,int*); 13ca44d042SBarry Smith EXTERN int MatGetRow_MPIBAIJ(Mat,int,int*,int**,Scalar**); 14ca44d042SBarry Smith EXTERN int MatRestoreRow_MPIBAIJ(Mat,int,int*,int**,Scalar**); 15d5b485abSSatish Balay 165615d1e5SSatish Balay #undef __FUNC__ 17b0a32e0cSBarry Smith #define __FUNC__ "MatCompressIndicesGeneral_MPIBAIJ" 18e757655aSSatish Balay static int MatCompressIndicesGeneral_MPIBAIJ(Mat C,int imax,IS *is_in,IS *is_out) 19df36731bSBarry Smith { 20df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)C->data; 21233c2ff6SSatish Balay int ierr,isz,bs = baij->bs,n,i,j,*idx,ival; 22233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 23233c2ff6SSatish Balay PetscTable gid1_lid1; 24233c2ff6SSatish Balay int tt, gid1, *nidx; 25233c2ff6SSatish Balay PetscTablePosition tpos; 26233c2ff6SSatish Balay #else 27233c2ff6SSatish Balay int Nbs,*nidx; 28f1af5d2fSBarry Smith PetscBT table; 29233c2ff6SSatish Balay #endif 30df36731bSBarry Smith 313a40ed3dSBarry Smith PetscFunctionBegin; 32233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 33233c2ff6SSatish Balay ierr = PetscTableCreate(baij->mbs,&gid1_lid1);CHKERRQ(ierr); 34233c2ff6SSatish Balay #else 35df36731bSBarry Smith Nbs = baij->Nbs; 36b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(int),&nidx);CHKERRQ(ierr); 376831982aSBarry Smith ierr = PetscBTCreate(Nbs,table);CHKERRQ(ierr); 38233c2ff6SSatish Balay #endif 39df36731bSBarry Smith for (i=0; i<imax; i++) { 40df36731bSBarry Smith isz = 0; 41233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 42233c2ff6SSatish Balay ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 43233c2ff6SSatish Balay #else 446831982aSBarry Smith ierr = PetscBTMemzero(Nbs,table);CHKERRQ(ierr); 45233c2ff6SSatish Balay #endif 4636f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 47b9b97703SBarry Smith ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 48e757655aSSatish Balay for (j=0; j<n ; j++) { 49df36731bSBarry Smith ival = idx[j]/bs; /* convert the indices into block indices */ 50233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 51233c2ff6SSatish Balay ierr = PetscTableFind(gid1_lid1,ival+1,&tt);CHKERRQ(ierr); 52233c2ff6SSatish Balay if (!tt) { 53233c2ff6SSatish Balay ierr = PetscTableAdd(gid1_lid1,ival+1,isz+1);CHKERRQ(ierr); 54233c2ff6SSatish Balay isz++; 55233c2ff6SSatish Balay } 56233c2ff6SSatish Balay #else 5729bbc08cSBarry Smith if (ival>Nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim"); 58f1af5d2fSBarry Smith if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;} 59233c2ff6SSatish Balay #endif 60df36731bSBarry Smith } 6136f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr); 62233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 63b0a32e0cSBarry Smith ierr = PetscMalloc((isz+1)*sizeof(int),&nidx);CHKERRQ(ierr); 64233c2ff6SSatish Balay ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 65233c2ff6SSatish Balay j = 0; 66233c2ff6SSatish Balay while (tpos) { 67233c2ff6SSatish Balay ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid1,&tt);CHKERRQ(ierr); 6829bbc08cSBarry Smith if (tt-- > isz) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than array-dim"); } 69233c2ff6SSatish Balay nidx[tt] = gid1 - 1; 70233c2ff6SSatish Balay j++; 71df36731bSBarry Smith } 7229bbc08cSBarry Smith if (j != isz) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"table error: jj != isz"); } 73233c2ff6SSatish Balay ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is_out+i));CHKERRQ(ierr); 74233c2ff6SSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 75233c2ff6SSatish Balay #else 76233c2ff6SSatish Balay ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is_out+i));CHKERRQ(ierr); 77233c2ff6SSatish Balay #endif 78233c2ff6SSatish Balay } 79233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 80233c2ff6SSatish Balay ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 81233c2ff6SSatish Balay #else 826831982aSBarry Smith ierr = PetscBTDestroy(table);CHKERRQ(ierr); 83606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 84233c2ff6SSatish Balay #endif 853a40ed3dSBarry Smith PetscFunctionReturn(0); 86df36731bSBarry Smith } 87df36731bSBarry Smith 885615d1e5SSatish Balay #undef __FUNC__ 89b0a32e0cSBarry Smith #define __FUNC__ "MatCompressIndicesSorted_MPIBAIJ" 90e757655aSSatish Balay static int MatCompressIndicesSorted_MPIBAIJ(Mat C,int imax,IS *is_in,IS *is_out) 91e757655aSSatish Balay { 92e757655aSSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)C->data; 93233c2ff6SSatish Balay int ierr,bs=baij->bs,i,j,k,val,n,*idx,*nidx,*idx_local; 94e757655aSSatish Balay PetscTruth flg; 95233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 96233c2ff6SSatish Balay int maxsz; 97233c2ff6SSatish Balay #else 98233c2ff6SSatish Balay int Nbs=baij->Nbs; 99233c2ff6SSatish Balay #endif 1003a40ed3dSBarry Smith PetscFunctionBegin; 101e757655aSSatish Balay for (i=0; i<imax; i++) { 102e757655aSSatish Balay ierr = ISSorted(is_in[i],&flg);CHKERRQ(ierr); 10329bbc08cSBarry Smith if (!flg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Indices are not sorted"); 104e757655aSSatish Balay } 105233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 106233c2ff6SSatish Balay /* Now check max size */ 107233c2ff6SSatish Balay for (i=0,maxsz=0; i<imax; i++) { 108233c2ff6SSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 109233c2ff6SSatish Balay ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 110ac355199SBarry Smith if (n%bs !=0) SETERRQ(1,"Indices are not block ordered"); 111233c2ff6SSatish Balay n = n/bs; /* The reduced index size */ 112233c2ff6SSatish Balay if (n > maxsz) maxsz = n; 113233c2ff6SSatish Balay } 114b0a32e0cSBarry Smith ierr = PetscMalloc((maxsz+1)*sizeof(int),&nidx);CHKERRQ(ierr); 115233c2ff6SSatish Balay #else 116b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(int),&nidx);CHKERRQ(ierr); 117233c2ff6SSatish Balay #endif 1183a40ed3dSBarry Smith /* Now check if the indices are in block order */ 119e757655aSSatish Balay for (i=0; i<imax; i++) { 120e757655aSSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 121b9b97703SBarry Smith ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 122ac355199SBarry Smith if (n%bs !=0) SETERRQ(1,"Indices are not block ordered"); 123e757655aSSatish Balay 124e757655aSSatish Balay n = n/bs; /* The reduced index size */ 125e757655aSSatish Balay idx_local = idx; 126e757655aSSatish Balay for (j=0; j<n ; j++) { 127e757655aSSatish Balay val = idx_local[0]; 128ac355199SBarry Smith if (val%bs != 0) SETERRQ(1,"Indices are not block ordered"); 129e757655aSSatish Balay for (k=0; k<bs; k++) { 130ac355199SBarry Smith if (val+k != idx_local[k]) SETERRQ(1,"Indices are not block ordered"); 131e757655aSSatish Balay } 132e757655aSSatish Balay nidx[j] = val/bs; 133e757655aSSatish Balay idx_local +=bs; 134e757655aSSatish Balay } 135e757655aSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr); 136029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,nidx,(is_out+i));CHKERRQ(ierr); 137e757655aSSatish Balay } 138606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 139233c2ff6SSatish Balay 1403a40ed3dSBarry Smith PetscFunctionReturn(0); 141e757655aSSatish Balay } 142e757655aSSatish Balay 1435615d1e5SSatish Balay #undef __FUNC__ 144b0a32e0cSBarry Smith #define __FUNC__ "MatExpandIndices_MPIBAIJ" 14536f4e84dSSatish Balay static int MatExpandIndices_MPIBAIJ(Mat C,int imax,IS *is_in,IS *is_out) 146df36731bSBarry Smith { 147df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)C->data; 148233c2ff6SSatish Balay int ierr,bs = baij->bs,n,i,j,k,*idx,*nidx; 149233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 150233c2ff6SSatish Balay int maxsz; 151233c2ff6SSatish Balay #else 152233c2ff6SSatish Balay int Nbs = baij->Nbs; 153233c2ff6SSatish Balay #endif 1543a40ed3dSBarry Smith PetscFunctionBegin; 155df36731bSBarry Smith 156233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 157233c2ff6SSatish Balay /* Now check max size */ 158233c2ff6SSatish Balay for (i=0,maxsz=0; i<imax; i++) { 159233c2ff6SSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 160233c2ff6SSatish Balay ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 161233c2ff6SSatish Balay if (n*bs > maxsz) maxsz = n*bs; 162233c2ff6SSatish Balay } 163b0a32e0cSBarry Smith ierr = PetscMalloc((maxsz+1)*sizeof(int),&nidx);CHKERRQ(ierr); 164233c2ff6SSatish Balay #else 165b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs*bs+1)*sizeof(int),&nidx);CHKERRQ(ierr); 166233c2ff6SSatish Balay #endif 167df36731bSBarry Smith 168df36731bSBarry Smith for (i=0; i<imax; i++) { 16936f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 170b9b97703SBarry Smith ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 171df36731bSBarry Smith for (j=0; j<n ; ++j){ 172df36731bSBarry Smith for (k=0; k<bs; k++) 173df36731bSBarry Smith nidx[j*bs+k] = idx[j]*bs+k; 174df36731bSBarry Smith } 17536f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr); 176329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n*bs,nidx,is_out+i);CHKERRQ(ierr); 177df36731bSBarry Smith } 178606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 1793a40ed3dSBarry Smith PetscFunctionReturn(0); 180df36731bSBarry Smith } 181df36731bSBarry Smith 182df36731bSBarry Smith 1835615d1e5SSatish Balay #undef __FUNC__ 184b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ" 185df36731bSBarry Smith int MatIncreaseOverlap_MPIBAIJ(Mat C,int imax,IS *is,int ov) 186d5b485abSSatish Balay { 187d5b485abSSatish Balay int i,ierr; 18836f4e84dSSatish Balay IS *is_new; 18936f4e84dSSatish Balay 1903a40ed3dSBarry Smith PetscFunctionBegin; 19182502324SSatish Balay ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr); 192df36731bSBarry Smith /* Convert the indices into block format */ 193e757655aSSatish Balay ierr = MatCompressIndicesGeneral_MPIBAIJ(C,imax,is,is_new);CHKERRQ(ierr); 19429bbc08cSBarry Smith if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");} 195d5b485abSSatish Balay for (i=0; i<ov; ++i) { 19636f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 197d5b485abSSatish Balay } 198ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);} 19936f4e84dSSatish Balay ierr = MatExpandIndices_MPIBAIJ(C,imax,is_new,is);CHKERRQ(ierr); 200ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);} 201606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 203d5b485abSSatish Balay } 204d5b485abSSatish Balay 205d5b485abSSatish Balay /* 206d5b485abSSatish Balay Sample message format: 207d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 208d5b485abSSatish Balay to index sets 1s[1], is[5] 209d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 210d5b485abSSatish Balay ----------- 211d5b485abSSatish Balay mesg [1] = 1 => is[1] 212d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 213d5b485abSSatish Balay ----------- 214d5b485abSSatish Balay mesg [5] = 5 => is[5] 215d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 216d5b485abSSatish Balay ----------- 217d5b485abSSatish Balay mesg [7] 218d5b485abSSatish Balay mesg [n] datas[1] 219d5b485abSSatish Balay ----------- 220d5b485abSSatish Balay mesg[n+1] 221d5b485abSSatish Balay mesg[m] data(is[5]) 222d5b485abSSatish Balay ----------- 223d5b485abSSatish Balay 224d5b485abSSatish Balay Notes: 225d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 226d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 227d5b485abSSatish Balay */ 2285615d1e5SSatish Balay #undef __FUNC__ 229b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Once" 230df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat C,int imax,IS *is) 231d5b485abSSatish Balay { 232df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 233d5b485abSSatish Balay int **idx,*n,*w1,*w2,*w3,*w4,*rtable,**data,len,*idx_i; 234df36731bSBarry Smith int size,rank,Mbs,i,j,k,ierr,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 235c7dd2924SSatish Balay int *ctr,*pa,*tmp,nrqr,*isz,*isz1,**xdata,**rbuf2; 236c7dd2924SSatish Balay int *onodes1,*olengths1,tag1,tag2,*onodes2,*olengths2; 237f1af5d2fSBarry Smith PetscBT *table; 238d5b485abSSatish Balay MPI_Comm comm; 239d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 240d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 241d5b485abSSatish Balay 2423a40ed3dSBarry Smith PetscFunctionBegin; 243d5b485abSSatish Balay comm = C->comm; 244d5b485abSSatish Balay size = c->size; 245d5b485abSSatish Balay rank = c->rank; 246df36731bSBarry Smith Mbs = c->Mbs; 247d5b485abSSatish Balay 248c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 249c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 250c7dd2924SSatish Balay 251df36731bSBarry Smith len = (imax+1)*sizeof(int*)+ (imax + Mbs)*sizeof(int); 25282502324SSatish Balay ierr = PetscMalloc(len,&idx);CHKERRQ(ierr); 253d5b485abSSatish Balay n = (int*)(idx + imax); 254d5b485abSSatish Balay rtable = n + imax; 255d5b485abSSatish Balay 256d5b485abSSatish Balay for (i=0; i<imax; i++) { 257d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 258b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 259d5b485abSSatish Balay } 260d5b485abSSatish Balay 261d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 262d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 263d5b485abSSatish Balay len = c->rowners[i+1]; 264d5b485abSSatish Balay for (; j<len; j++) { 265d5b485abSSatish Balay rtable[j] = i; 266d5b485abSSatish Balay } 267d5b485abSSatish Balay } 268d5b485abSSatish Balay 269d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 270d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 27182502324SSatish Balay ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr);/* mesg size */ 272d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 273d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 274d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 275549d3d68SSatish Balay ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 276d5b485abSSatish Balay for (i=0; i<imax; i++) { 277549d3d68SSatish Balay ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 278d5b485abSSatish Balay idx_i = idx[i]; 279d5b485abSSatish Balay len = n[i]; 280d5b485abSSatish Balay for (j=0; j<len; j++) { 281d5b485abSSatish Balay row = idx_i[j]; 2826b41c64dSBarry Smith if (row < 0) { 28329bbc08cSBarry Smith SETERRQ(1,"Index set cannot have negative entries"); 2846b41c64dSBarry Smith } 285d5b485abSSatish Balay proc = rtable[row]; 286d5b485abSSatish Balay w4[proc]++; 287d5b485abSSatish Balay } 288d5b485abSSatish Balay for (j=0; j<size; j++){ 289d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 290d5b485abSSatish Balay } 291d5b485abSSatish Balay } 292d5b485abSSatish Balay 293d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 294d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 295d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 296d5b485abSSatish Balay w3[rank] = 0; 297d5b485abSSatish Balay for (i=0; i<size; i++) { 298d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 299d5b485abSSatish Balay } 300d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 301b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); 302d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 303d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 304d5b485abSSatish Balay } 305d5b485abSSatish Balay 306d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 307d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 308d5b485abSSatish Balay j = pa[i]; 309d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 310d5b485abSSatish Balay msz += w1[j]; 311d5b485abSSatish Balay } 312d5b485abSSatish Balay 313c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 314c7dd2924SSatish Balay ierr = PetscGatherNumberOfMessages(comm,nrqs,w2,w1,&nrqr);CHKERRQ(ierr); 315c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 316d5b485abSSatish Balay 317c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 318c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 319d5b485abSSatish Balay 320d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 321d5b485abSSatish Balay len = 2*size*sizeof(int*) + (size+msz)*sizeof(int); 32282502324SSatish Balay ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr); 323d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 324549d3d68SSatish Balay ierr = PetscMemzero(outdat,2*size*sizeof(int*));CHKERRQ(ierr); 325d5b485abSSatish Balay tmp = (int*)(outdat + 2*size); 326d5b485abSSatish Balay ctr = tmp + msz; 327d5b485abSSatish Balay 328d5b485abSSatish Balay { 329d5b485abSSatish Balay int *iptr = tmp,ict = 0; 330d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 331d5b485abSSatish Balay j = pa[i]; 332d5b485abSSatish Balay iptr += ict; 333d5b485abSSatish Balay outdat[j] = iptr; 334d5b485abSSatish Balay ict = w1[j]; 335d5b485abSSatish Balay } 336d5b485abSSatish Balay } 337d5b485abSSatish Balay 338d5b485abSSatish Balay /* Form the outgoing messages */ 339d5b485abSSatish Balay /*plug in the headers*/ 340d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 341d5b485abSSatish Balay j = pa[i]; 342d5b485abSSatish Balay outdat[j][0] = 0; 343549d3d68SSatish Balay ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 344d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 345d5b485abSSatish Balay } 346d5b485abSSatish Balay 347d5b485abSSatish Balay /* Memory for doing local proc's work*/ 348d5b485abSSatish Balay { 349d5b485abSSatish Balay int *d_p; 350d5b485abSSatish Balay char *t_p; 351d5b485abSSatish Balay 352f1af5d2fSBarry Smith len = (imax)*(sizeof(PetscBT) + sizeof(int*)+ sizeof(int)) + 3532a8f2162SSatish Balay (Mbs)*imax*sizeof(int) + (Mbs/BITSPERBYTE+1)*imax*sizeof(char) + 1; 35482502324SSatish Balay ierr = PetscMalloc(len,&table);CHKERRQ(ierr); 355549d3d68SSatish Balay ierr = PetscMemzero(table,len);CHKERRQ(ierr); 356d5b485abSSatish Balay data = (int **)(table + imax); 357bbd702dbSSatish Balay isz = (int *)(data + imax); 358bbd702dbSSatish Balay d_p = (int *)(isz + imax); 359bbd702dbSSatish Balay t_p = (char *)(d_p + Mbs*imax); 360bbd702dbSSatish Balay for (i=0; i<imax; i++) { 3612a8f2162SSatish Balay table[i] = t_p + (Mbs/BITSPERBYTE+1)*i; 362bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 363d5b485abSSatish Balay } 364d5b485abSSatish Balay } 365d5b485abSSatish Balay 366d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 367d5b485abSSatish Balay { 368d5b485abSSatish Balay int n_i,*data_i,isz_i,*outdat_j,ctr_j; 369f1af5d2fSBarry Smith PetscBT table_i; 370d5b485abSSatish Balay 371d5b485abSSatish Balay for (i=0; i<imax; i++) { 372549d3d68SSatish Balay ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 373d5b485abSSatish Balay n_i = n[i]; 374d5b485abSSatish Balay table_i = table[i]; 375d5b485abSSatish Balay idx_i = idx[i]; 376d5b485abSSatish Balay data_i = data[i]; 377d5b485abSSatish Balay isz_i = isz[i]; 378d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 379d5b485abSSatish Balay row = idx_i[j]; 380d5b485abSSatish Balay proc = rtable[row]; 381d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 382d5b485abSSatish Balay ctr[proc]++; 383d5b485abSSatish Balay *ptr[proc] = row; 384d5b485abSSatish Balay ptr[proc]++; 385d5b485abSSatish Balay } 386d5b485abSSatish Balay else { /* Update the local table */ 387f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 388d5b485abSSatish Balay } 389d5b485abSSatish Balay } 390d5b485abSSatish Balay /* Update the headers for the current IS */ 391d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 392d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 393d5b485abSSatish Balay outdat_j = outdat[j]; 394d5b485abSSatish Balay k = ++outdat_j[0]; 395d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 396d5b485abSSatish Balay outdat_j[2*k-1] = i; 397d5b485abSSatish Balay } 398d5b485abSSatish Balay } 399d5b485abSSatish Balay isz[i] = isz_i; 400d5b485abSSatish Balay } 401d5b485abSSatish Balay } 402d5b485abSSatish Balay 403d5b485abSSatish Balay /* Now post the sends */ 404b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 405d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 406d5b485abSSatish Balay j = pa[i]; 407c7dd2924SSatish Balay ierr = MPI_Isend(outdat[j],w1[j],MPI_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 408d5b485abSSatish Balay } 409d5b485abSSatish Balay 410d5b485abSSatish Balay /* No longer need the original indices*/ 411d5b485abSSatish Balay for (i=0; i<imax; ++i) { 412d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 413d5b485abSSatish Balay } 414606d414cSSatish Balay ierr = PetscFree(idx);CHKERRQ(ierr); 415d5b485abSSatish Balay 416d5b485abSSatish Balay for (i=0; i<imax; ++i) { 417d5b485abSSatish Balay ierr = ISDestroy(is[i]);CHKERRQ(ierr); 418d5b485abSSatish Balay } 419d5b485abSSatish Balay 420d5b485abSSatish Balay /* Do Local work*/ 421df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 422d5b485abSSatish Balay 423d5b485abSSatish Balay /* Receive messages*/ 424b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 425c7dd2924SSatish Balay ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr); 426d5b485abSSatish Balay 427b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 428ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr); 429d5b485abSSatish Balay 430d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 431606d414cSSatish Balay ierr = PetscFree(outdat);CHKERRQ(ierr); 432606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 433d5b485abSSatish Balay 43482502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(int *),&xdata);CHKERRQ(ierr); 435b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(int),&isz1);CHKERRQ(ierr); 436df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 437606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 438d5b485abSSatish Balay 439d5b485abSSatish Balay /* Send the data back*/ 440d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 441d5b485abSSatish Balay { 442c7dd2924SSatish Balay int *rw1; 443d5b485abSSatish Balay 444c7dd2924SSatish Balay ierr = PetscMalloc(size*sizeof(int),&rw1);CHKERRQ(ierr); 445c7dd2924SSatish Balay ierr = PetscMemzero(rw1,size*sizeof(int));CHKERRQ(ierr); 446c7dd2924SSatish Balay 447d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 448d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 449c7dd2924SSatish Balay if (proc != onodes1[i]) SETERRQ(1,"MPI_SOURCE mismatch"); 450d5b485abSSatish Balay rw1[proc] = isz1[i]; 451d5b485abSSatish Balay } 452d5b485abSSatish Balay 453c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 454c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 455d5b485abSSatish Balay 456c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 457c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 458c7dd2924SSatish Balay PetscFree(rw1); 459c7dd2924SSatish Balay } 460c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 461c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 462d5b485abSSatish Balay 463d5b485abSSatish Balay /* Now post the sends */ 464b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 465d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 466d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 467c7dd2924SSatish Balay ierr = MPI_Isend(xdata[i],isz1[i],MPI_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 468d5b485abSSatish Balay } 469d5b485abSSatish Balay 470d5b485abSSatish Balay /* receive work done on other processors*/ 471d5b485abSSatish Balay { 472d5b485abSSatish Balay int index,is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 473f1af5d2fSBarry Smith PetscBT table_i; 474d5b485abSSatish Balay MPI_Status *status2; 475d5b485abSSatish Balay 476b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 477d5b485abSSatish Balay 478d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 479ca161407SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&index,status2+i);CHKERRQ(ierr); 480d5b485abSSatish Balay /* Process the message*/ 481d5b485abSSatish Balay rbuf2_i = rbuf2[index]; 482d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 483d5b485abSSatish Balay jmax = rbuf2[index][0]; 484d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 485d5b485abSSatish Balay max = rbuf2_i[2*j]; 486d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 487d5b485abSSatish Balay isz_i = isz[is_no]; 488d5b485abSSatish Balay data_i = data[is_no]; 489d5b485abSSatish Balay table_i = table[is_no]; 490d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 491d5b485abSSatish Balay row = rbuf2_i[ct1]; 492f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 493d5b485abSSatish Balay } 494d5b485abSSatish Balay isz[is_no] = isz_i; 495d5b485abSSatish Balay } 496d5b485abSSatish Balay } 497ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr); 498606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 499d5b485abSSatish Balay } 500d5b485abSSatish Balay 501d5b485abSSatish Balay for (i=0; i<imax; ++i) { 502029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr); 503d5b485abSSatish Balay } 504d5b485abSSatish Balay 505c7dd2924SSatish Balay 506c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 507c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 508c7dd2924SSatish Balay 509606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 510606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 511606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 512606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 513606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 514606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 515606d414cSSatish Balay ierr = PetscFree(table);CHKERRQ(ierr); 516606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 517606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 518606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 519606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 520606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 5213a40ed3dSBarry Smith PetscFunctionReturn(0); 522d5b485abSSatish Balay } 523d5b485abSSatish Balay 5245615d1e5SSatish Balay #undef __FUNC__ 525b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Local" 526d5b485abSSatish Balay /* 527df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 528d5b485abSSatish Balay the work on the local processor. 529d5b485abSSatish Balay 530d5b485abSSatish Balay Inputs: 531df36731bSBarry Smith C - MAT_MPIBAIJ; 532d5b485abSSatish Balay imax - total no of index sets processed at a time; 533df36731bSBarry Smith table - an array of char - size = Mbs bits. 534d5b485abSSatish Balay 535d5b485abSSatish Balay Output: 536d5b485abSSatish Balay isz - array containing the count of the solution elements correspondign 537d5b485abSSatish Balay to each index set; 538d5b485abSSatish Balay data - pointer to the solutions 539d5b485abSSatish Balay */ 540f1af5d2fSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,PetscBT *table,int *isz,int **data) 541d5b485abSSatish Balay { 542df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 543d5b485abSSatish Balay Mat A = c->A,B = c->B; 544df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 545df36731bSBarry Smith int start,end,val,max,rstart,cstart,*ai,*aj; 546d5b485abSSatish Balay int *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 547f1af5d2fSBarry Smith PetscBT table_i; 548d5b485abSSatish Balay 5493a40ed3dSBarry Smith PetscFunctionBegin; 550d5b485abSSatish Balay rstart = c->rstart; 551d5b485abSSatish Balay cstart = c->cstart; 552d5b485abSSatish Balay ai = a->i; 553df36731bSBarry Smith aj = a->j; 554d5b485abSSatish Balay bi = b->i; 555df36731bSBarry Smith bj = b->j; 556d5b485abSSatish Balay garray = c->garray; 557d5b485abSSatish Balay 558d5b485abSSatish Balay 559d5b485abSSatish Balay for (i=0; i<imax; i++) { 560d5b485abSSatish Balay data_i = data[i]; 561d5b485abSSatish Balay table_i = table[i]; 562d5b485abSSatish Balay isz_i = isz[i]; 563d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 564d5b485abSSatish Balay row = data_i[j] - rstart; 565d5b485abSSatish Balay start = ai[row]; 566d5b485abSSatish Balay end = ai[row+1]; 567d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 568df36731bSBarry Smith val = aj[k] + cstart; 569f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 570d5b485abSSatish Balay } 571d5b485abSSatish Balay start = bi[row]; 572d5b485abSSatish Balay end = bi[row+1]; 573d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 574df36731bSBarry Smith val = garray[bj[k]]; 575f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 576d5b485abSSatish Balay } 577d5b485abSSatish Balay } 578d5b485abSSatish Balay isz[i] = isz_i; 579d5b485abSSatish Balay } 5803a40ed3dSBarry Smith PetscFunctionReturn(0); 581d5b485abSSatish Balay } 5825615d1e5SSatish Balay #undef __FUNC__ 583b0a32e0cSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Receive" 584d5b485abSSatish Balay /* 585df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 586d5b485abSSatish Balay and return the output 587d5b485abSSatish Balay 588d5b485abSSatish Balay Input: 589d5b485abSSatish Balay C - the matrix 590d5b485abSSatish Balay nrqr - no of messages being processed. 591d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 592d5b485abSSatish Balay 593d5b485abSSatish Balay Output: 594d5b485abSSatish Balay xdata - array of messages to be sent back 595d5b485abSSatish Balay isz1 - size of each message 596d5b485abSSatish Balay 597d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 598d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 599d5b485abSSatish Balay rather then all previous rows as it is now where a single large chunck of 600d5b485abSSatish Balay memory is used. 601d5b485abSSatish Balay 602d5b485abSSatish Balay */ 603701b8814SBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf,int **xdata,int * isz1) 604d5b485abSSatish Balay { 605df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 606d5b485abSSatish Balay Mat A = c->A,B = c->B; 607df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 608df36731bSBarry Smith int rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 609d5b485abSSatish Balay int row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 610df36731bSBarry Smith int val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 611bbd702dbSSatish Balay int *rbuf_i,kmax,rbuf_0,ierr; 612f1af5d2fSBarry Smith PetscBT xtable; 613d5b485abSSatish Balay 6143a40ed3dSBarry Smith PetscFunctionBegin; 615d5b485abSSatish Balay rank = c->rank; 616df36731bSBarry Smith Mbs = c->Mbs; 617d5b485abSSatish Balay rstart = c->rstart; 618d5b485abSSatish Balay cstart = c->cstart; 619d5b485abSSatish Balay ai = a->i; 620df36731bSBarry Smith aj = a->j; 621d5b485abSSatish Balay bi = b->i; 622df36731bSBarry Smith bj = b->j; 623d5b485abSSatish Balay garray = c->garray; 624d5b485abSSatish Balay 625d5b485abSSatish Balay 626d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 627d5b485abSSatish Balay rbuf_i = rbuf[i]; 628d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 629d5b485abSSatish Balay ct += rbuf_0; 630d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 631d5b485abSSatish Balay } 632d5b485abSSatish Balay 633701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 634701b8814SBarry Smith else max1 = 1; 635d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 63682502324SSatish Balay ierr = PetscMalloc(mem_estimate*sizeof(int),&xdata[0]);CHKERRQ(ierr); 637d5b485abSSatish Balay ++no_malloc; 6386831982aSBarry Smith ierr = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr); 639549d3d68SSatish Balay ierr = PetscMemzero(isz1,nrqr*sizeof(int));CHKERRQ(ierr); 640d5b485abSSatish Balay 641d5b485abSSatish Balay ct3 = 0; 642d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 643d5b485abSSatish Balay rbuf_i = rbuf[i]; 644d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 645d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 646d5b485abSSatish Balay ct2 = ct1; 647d5b485abSSatish Balay ct3 += ct1; 648d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 6496831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 650d5b485abSSatish Balay oct2 = ct2; 651d5b485abSSatish Balay kmax = rbuf_i[2*j]; 652d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 653d5b485abSSatish Balay row = rbuf_i[ct1]; 654f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 655d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 656d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 65782502324SSatish Balay ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr); 658549d3d68SSatish Balay ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 659606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 660d5b485abSSatish Balay xdata[0] = tmp; 661d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 662d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 663d5b485abSSatish Balay } 664d5b485abSSatish Balay xdata[i][ct2++] = row; 665d5b485abSSatish Balay ct3++; 666d5b485abSSatish Balay } 667d5b485abSSatish Balay } 668d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 669d5b485abSSatish Balay row = xdata[i][k] - rstart; 670d5b485abSSatish Balay start = ai[row]; 671d5b485abSSatish Balay end = ai[row+1]; 672d5b485abSSatish Balay for (l=start; l<end; l++) { 673df36731bSBarry Smith val = aj[l] + cstart; 674f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 675d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 676d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 67782502324SSatish Balay ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr); 678549d3d68SSatish Balay ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 679606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 680d5b485abSSatish Balay xdata[0] = tmp; 681d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 682d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 683d5b485abSSatish Balay } 684d5b485abSSatish Balay xdata[i][ct2++] = val; 685d5b485abSSatish Balay ct3++; 686d5b485abSSatish Balay } 687d5b485abSSatish Balay } 688d5b485abSSatish Balay start = bi[row]; 689d5b485abSSatish Balay end = bi[row+1]; 690d5b485abSSatish Balay for (l=start; l<end; l++) { 691df36731bSBarry Smith val = garray[bj[l]]; 692f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 693d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 694d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 69582502324SSatish Balay ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr); 696549d3d68SSatish Balay ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 697606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 698d5b485abSSatish Balay xdata[0] = tmp; 699d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 700d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 701d5b485abSSatish Balay } 702d5b485abSSatish Balay xdata[i][ct2++] = val; 703d5b485abSSatish Balay ct3++; 704d5b485abSSatish Balay } 705d5b485abSSatish Balay } 706d5b485abSSatish Balay } 707d5b485abSSatish Balay /* Update the header*/ 708d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 709d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 710d5b485abSSatish Balay } 711d5b485abSSatish Balay xdata[i][0] = rbuf_0; 712d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 713d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 714d5b485abSSatish Balay } 7156831982aSBarry Smith ierr = PetscBTDestroy(xtable);CHKERRQ(ierr); 716b0a32e0cSBarry Smith PetscLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d, no of mallocs = %d\n",rank,mem_estimate,ct3,no_malloc); 7173a40ed3dSBarry Smith PetscFunctionReturn(0); 718d5b485abSSatish Balay } 719d5b485abSSatish Balay 7207b2a1423SBarry Smith static int MatGetSubMatrices_MPIBAIJ_local(Mat,int,IS *,IS *,MatReuse,Mat *); 721a2feddc5SSatish Balay 7225615d1e5SSatish Balay #undef __FUNC__ 723b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrices_MPIBAIJ" 724f1af5d2fSBarry Smith int MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat **submat) 725d5b485abSSatish Balay { 72636f4e84dSSatish Balay IS *isrow_new,*iscol_new; 727cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 728cf4f527aSSatish Balay int nmax,nstages_local,nstages,i,pos,max_no,ierr; 729a2feddc5SSatish Balay 7303a40ed3dSBarry Smith PetscFunctionBegin; 7313a40ed3dSBarry Smith /* The compression and expansion should be avoided. Does'nt point 732a2feddc5SSatish Balay out errors might change the indices hence buggey */ 733a2feddc5SSatish Balay 734b0a32e0cSBarry Smith ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&isrow_new);CHKERRQ(ierr); 73536f4e84dSSatish Balay iscol_new = isrow_new + ismax; 736e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C,ismax,isrow,isrow_new);CHKERRQ(ierr); 737e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C,ismax,iscol,iscol_new);CHKERRQ(ierr); 738cf4f527aSSatish Balay 739cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 740cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 74182502324SSatish Balay ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 742cf4f527aSSatish Balay } 743cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 744cf4f527aSSatish Balay nmax = 20*1000000 / (c->Nbs * sizeof(int)); 745329f5518SBarry Smith if (!nmax) nmax = 1; 746cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 747cf4f527aSSatish Balay 748cf4f527aSSatish Balay /* Make sure every porcessor loops through the nstages */ 749ca161407SBarry Smith ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 750cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 751cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 752cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 753cf4f527aSSatish Balay else max_no = ismax-pos; 754cf4f527aSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr); 755cf4f527aSSatish Balay pos += max_no; 756cf4f527aSSatish Balay } 75736f4e84dSSatish Balay 75836f4e84dSSatish Balay for (i=0; i<ismax; i++) { 759ca161407SBarry Smith ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr); 760ca161407SBarry Smith ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr); 76136f4e84dSSatish Balay } 762606d414cSSatish Balay ierr = PetscFree(isrow_new);CHKERRQ(ierr); 7633a40ed3dSBarry Smith PetscFunctionReturn(0); 764a2feddc5SSatish Balay } 765a2feddc5SSatish Balay 766233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 767233c2ff6SSatish Balay #undef __FUNC__ 768233c2ff6SSatish Balay #define __FUNC__ "PetscGetProc" 769233c2ff6SSatish Balay int PetscGetProc(const int gid, const int numprocs, const int proc_gnode[], int *proc) 770233c2ff6SSatish Balay { 771233c2ff6SSatish Balay int nGlobalNd = proc_gnode[numprocs]; 772233c2ff6SSatish Balay int fproc = (int) ((float)gid * (float)numprocs / (float)nGlobalNd + 0.5); 773233c2ff6SSatish Balay 774233c2ff6SSatish Balay PetscFunctionBegin; 77529bbc08cSBarry Smith /* if(fproc < 0) SETERRQ(1,"fproc < 0");*/ 776233c2ff6SSatish Balay if (fproc > numprocs) fproc = numprocs; 777233c2ff6SSatish Balay while (gid < proc_gnode[fproc] || gid >= proc_gnode[fproc+1]) { 778233c2ff6SSatish Balay if (gid < proc_gnode[fproc]) fproc--; 779233c2ff6SSatish Balay else fproc++; 780233c2ff6SSatish Balay } 78129bbc08cSBarry Smith /* if(fproc<0 || fproc>=numprocs) { SETERRQ(1,"fproc < 0 || fproc >= numprocs"); }*/ 782233c2ff6SSatish Balay *proc = fproc; 783233c2ff6SSatish Balay PetscFunctionReturn(0); 784233c2ff6SSatish Balay } 785233c2ff6SSatish Balay #endif 786233c2ff6SSatish Balay 787a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 7885615d1e5SSatish Balay #undef __FUNC__ 789b0a32e0cSBarry Smith #define __FUNC__ "MatGetSubMatrices_MPIBAIJ_local" 7903eda8832SBarry Smith static int MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats) 791a2feddc5SSatish Balay { 792df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 793cf4f527aSSatish Balay Mat A = c->A; 794df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 795233c2ff6SSatish Balay int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,start,end,size; 796233c2ff6SSatish Balay int **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,ierr,**rbuf1,row,proc; 797c7dd2924SSatish Balay int nrqs,msz,**ptr,index,*req_size,*ctr,*pa,*tmp,tcol,nrqr; 798233c2ff6SSatish Balay int **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 799233c2ff6SSatish Balay int **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 800233c2ff6SSatish Balay int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 801233c2ff6SSatish Balay int bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 80236f4e84dSSatish Balay int cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 80334e6ae18SSatish Balay int *bmap = c->garray,ctmp,rstart=c->rstart,tag0,tag1,tag2,tag3; 804d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 805d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 806d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 807d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 808d5b485abSSatish Balay MPI_Comm comm; 8093eda8832SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 8103eda8832SBarry Smith MatScalar *a_a=a->a,*b_a=b->a; 8110f5bd95cSBarry Smith PetscTruth flag; 812c7dd2924SSatish Balay int *onodes1,*olengths1; 813c7dd2924SSatish Balay 81480d608c0SSatish Balay #if defined (PETSC_USE_CTABLE) 815233c2ff6SSatish Balay int tt; 816233c2ff6SSatish Balay PetscTable *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1; 817233c2ff6SSatish Balay #else 818233c2ff6SSatish Balay int **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 819233c2ff6SSatish Balay #endif 820d5b485abSSatish Balay 8213a40ed3dSBarry Smith PetscFunctionBegin; 822d5b485abSSatish Balay comm = C->comm; 82334e6ae18SSatish Balay tag0 = C->tag; 824d5b485abSSatish Balay size = c->size; 825d5b485abSSatish Balay rank = c->rank; 826d5b485abSSatish Balay 82734e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 82834e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 82934e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 83034e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 83134e6ae18SSatish Balay 832d5b485abSSatish Balay /* Check if the col indices are sorted */ 833d5b485abSSatish Balay for (i=0; i<ismax; i++) { 834888f2ed8SSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr); 83529bbc08cSBarry Smith if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 836d5b485abSSatish Balay } 837d5b485abSSatish Balay 838233c2ff6SSatish Balay len = (2*ismax+1)*(sizeof(int*)+ sizeof(int)); 839233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 840233c2ff6SSatish Balay len += (Mbs+1)*sizeof(int); 841233c2ff6SSatish Balay #endif 84282502324SSatish Balay ierr = PetscMalloc(len,&irow);CHKERRQ(ierr); 843d5b485abSSatish Balay icol = irow + ismax; 844d5b485abSSatish Balay nrow = (int*)(icol + ismax); 845d5b485abSSatish Balay ncol = nrow + ismax; 846233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 847d5b485abSSatish Balay rtable = ncol + ismax; 848d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 849d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 850d5b485abSSatish Balay jmax = c->rowners[i+1]; 851d5b485abSSatish Balay for (; j<jmax; j++) { 852d5b485abSSatish Balay rtable[j] = i; 853d5b485abSSatish Balay } 854d5b485abSSatish Balay } 855233c2ff6SSatish Balay #endif 856233c2ff6SSatish Balay 857233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 858233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 859233c2ff6SSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 860233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 861233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 862233c2ff6SSatish Balay } 863d5b485abSSatish Balay 864d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 865d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 86682502324SSatish Balay ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr); /* mesg size */ 867d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 868d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 869d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 870549d3d68SSatish Balay ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 871d5b485abSSatish Balay for (i=0; i<ismax; i++) { 872549d3d68SSatish Balay ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 873d5b485abSSatish Balay jmax = nrow[i]; 874d5b485abSSatish Balay irow_i = irow[i]; 875d5b485abSSatish Balay for (j=0; j<jmax; j++) { 876d5b485abSSatish Balay row = irow_i[j]; 877233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 878233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 879233c2ff6SSatish Balay #else 880d5b485abSSatish Balay proc = rtable[row]; 881233c2ff6SSatish Balay #endif 882d5b485abSSatish Balay w4[proc]++; 883d5b485abSSatish Balay } 884d5b485abSSatish Balay for (j=0; j<size; j++) { 885d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 886d5b485abSSatish Balay } 887d5b485abSSatish Balay } 888d5b485abSSatish Balay 889d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 890e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 891d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 892d5b485abSSatish Balay w3[rank] = 0; 893d5b485abSSatish Balay for (i=0; i<size; i++) { 894d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 895d5b485abSSatish Balay } 896b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); /*(proc -array)*/ 897d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 898d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 899d5b485abSSatish Balay } 900d5b485abSSatish Balay 901d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 902d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 903d5b485abSSatish Balay j = pa[i]; 904d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 905d5b485abSSatish Balay msz += w1[j]; 906d5b485abSSatish Balay } 907d5b485abSSatish Balay 908c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 909c7dd2924SSatish Balay ierr = PetscGatherNumberOfMessages(comm,nrqs,w2,w1,&nrqr);CHKERRQ(ierr); 910c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 911d5b485abSSatish Balay 912c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 913c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 914c7dd2924SSatish Balay 915c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 916c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 917d5b485abSSatish Balay 918d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 919d5b485abSSatish Balay len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 92082502324SSatish Balay ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr); 921d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 922549d3d68SSatish Balay ierr = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr); 923d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 924d5b485abSSatish Balay tmp = (int*)(ptr + size); 925d5b485abSSatish Balay ctr = tmp + 2*msz; 926d5b485abSSatish Balay 927d5b485abSSatish Balay { 928d5b485abSSatish Balay int *iptr = tmp,ict = 0; 929d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 930d5b485abSSatish Balay j = pa[i]; 931d5b485abSSatish Balay iptr += ict; 932d5b485abSSatish Balay sbuf1[j] = iptr; 933d5b485abSSatish Balay ict = w1[j]; 934d5b485abSSatish Balay } 935d5b485abSSatish Balay } 936d5b485abSSatish Balay 937d5b485abSSatish Balay /* Form the outgoing messages */ 938d5b485abSSatish Balay /* Initialise the header space */ 939d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 940d5b485abSSatish Balay j = pa[i]; 941d5b485abSSatish Balay sbuf1[j][0] = 0; 942549d3d68SSatish Balay ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 943d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 944d5b485abSSatish Balay } 945d5b485abSSatish Balay 946d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 947d5b485abSSatish Balay for (i=0; i<ismax; i++) { 948549d3d68SSatish Balay ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 949d5b485abSSatish Balay irow_i = irow[i]; 950d5b485abSSatish Balay jmax = nrow[i]; 951d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 952d5b485abSSatish Balay row = irow_i[j]; 953233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 954233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 955233c2ff6SSatish Balay #else 956d5b485abSSatish Balay proc = rtable[row]; 957233c2ff6SSatish Balay #endif 958d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 959d5b485abSSatish Balay ctr[proc]++; 960d5b485abSSatish Balay *ptr[proc] = row; 961d5b485abSSatish Balay ptr[proc]++; 962d5b485abSSatish Balay } 963d5b485abSSatish Balay } 964d5b485abSSatish Balay /* Update the headers for the current IS */ 965d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 966*06ef35abSSatish Balay if ((ctr_j = ctr[j])) { 967d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 968d5b485abSSatish Balay k = ++sbuf1_j[0]; 969d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 970d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 971d5b485abSSatish Balay } 972d5b485abSSatish Balay } 973d5b485abSSatish Balay } 974d5b485abSSatish Balay 975d5b485abSSatish Balay /* Now post the sends */ 976b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 977d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 978d5b485abSSatish Balay j = pa[i]; 979ca161407SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 980d5b485abSSatish Balay } 981d5b485abSSatish Balay 982d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 983b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 98482502324SSatish Balay ierr = PetscMalloc((nrqs+1)*sizeof(int *),&rbuf2);CHKERRQ(ierr); 985d5b485abSSatish Balay rbuf2[0] = tmp + msz; 986d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 987d5b485abSSatish Balay j = pa[i]; 988d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 989d5b485abSSatish Balay } 990d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 991d5b485abSSatish Balay j = pa[i]; 992ca161407SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 993d5b485abSSatish Balay } 994d5b485abSSatish Balay 995d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 996d5b485abSSatish Balay 997d5b485abSSatish Balay /* Receive messages*/ 998b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 999b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 1000d5b485abSSatish Balay len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*); 100182502324SSatish Balay ierr = PetscMalloc(len,&sbuf2);CHKERRQ(ierr); 1002d5b485abSSatish Balay req_size = (int*)(sbuf2 + nrqr); 1003d5b485abSSatish Balay req_source = req_size + nrqr; 1004d5b485abSSatish Balay 1005d5b485abSSatish Balay { 1006df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 1007ca161407SBarry Smith int *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 1008d5b485abSSatish Balay 1009d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1010ca161407SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&index,r_status1+i);CHKERRQ(ierr); 1011d5b485abSSatish Balay req_size[index] = 0; 1012d5b485abSSatish Balay rbuf1_i = rbuf1[index]; 1013d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 1014ca161407SBarry Smith ierr = MPI_Get_count(r_status1+i,MPI_INT,&end);CHKERRQ(ierr); 101582502324SSatish Balay ierr = PetscMalloc(end*sizeof(int),&sbuf2[index]);CHKERRQ(ierr); 1016d5b485abSSatish Balay sbuf2_i = sbuf2[index]; 1017d5b485abSSatish Balay for (j=start; j<end; j++) { 1018d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 1019d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 1020d5b485abSSatish Balay sbuf2_i[j] = ncols; 1021d5b485abSSatish Balay req_size[index] += ncols; 1022d5b485abSSatish Balay } 1023d5b485abSSatish Balay req_source[index] = r_status1[i].MPI_SOURCE; 1024d5b485abSSatish Balay /* form the header */ 1025d5b485abSSatish Balay sbuf2_i[0] = req_size[index]; 1026d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 1027ca161407SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPI_INT,req_source[index],tag1,comm,s_waits2+i);CHKERRQ(ierr); 1028d5b485abSSatish Balay } 1029d5b485abSSatish Balay } 1030606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 1031606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 1032d5b485abSSatish Balay 1033d5b485abSSatish Balay /* recv buffer sizes */ 1034d5b485abSSatish Balay /* Receive messages*/ 1035d5b485abSSatish Balay 1036b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);CHKERRQ(ierr); 1037b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr); 1038b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 1039b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 1040b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 1041d5b485abSSatish Balay 1042d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1043ca161407SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&index,r_status2+i);CHKERRQ(ierr); 104482502324SSatish Balay ierr = PetscMalloc(rbuf2[index][0]*sizeof(int),&rbuf3[index]);CHKERRQ(ierr); 104582502324SSatish Balay ierr = PetscMalloc(rbuf2[index][0]*bs2*sizeof(MatScalar),&rbuf4[index]);CHKERRQ(ierr); 1046ca161407SBarry Smith ierr = MPI_Irecv(rbuf3[index],rbuf2[index][0],MPI_INT, 1047ca161407SBarry Smith r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+index);CHKERRQ(ierr); 10483eda8832SBarry Smith ierr = MPI_Irecv(rbuf4[index],rbuf2[index][0]*bs2,MPIU_MATSCALAR, 1049ca161407SBarry Smith r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+index);CHKERRQ(ierr); 1050d5b485abSSatish Balay } 1051606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 1052606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 1053d5b485abSSatish Balay 1054d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 1055b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 1056b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 1057d5b485abSSatish Balay 1058ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 1059ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 1060606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 1061606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 1062606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 1063606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 1064d5b485abSSatish Balay 1065d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 106682502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(int *),&sbuf_aj);CHKERRQ(ierr); 1067d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 106882502324SSatish Balay ierr = PetscMalloc((j+1)*sizeof(int),&sbuf_aj[0]);CHKERRQ(ierr); 1069d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1070d5b485abSSatish Balay 1071b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 1072d5b485abSSatish Balay { 1073d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1074d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1075d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 1076d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 1077d5b485abSSatish Balay ct2 = 0; 1078d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1079d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 1080d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1081e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1082d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1083d5b485abSSatish Balay ncols = nzA + nzB; 1084d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1085d5b485abSSatish Balay 1086d5b485abSSatish Balay /* load the column indices for this row into cols*/ 1087d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 1088d5b485abSSatish Balay for (l=0; l<nzB; l++) { 1089d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1090d5b485abSSatish Balay else break; 1091d5b485abSSatish Balay } 1092d5b485abSSatish Balay imark = l; 1093d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1094d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1095d5b485abSSatish Balay ct2 += ncols; 1096d5b485abSSatish Balay } 1097d5b485abSSatish Balay } 1098ca161407SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 1099d5b485abSSatish Balay } 1100d5b485abSSatish Balay } 1101b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 1102b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 1103d5b485abSSatish Balay 1104d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 110582502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr); 1106d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 110782502324SSatish Balay ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr); 1108a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1109d5b485abSSatish Balay 1110b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 1111d5b485abSSatish Balay { 1112d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1113d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1114d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 1115d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 1116d5b485abSSatish Balay ct2 = 0; 1117d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1118d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 1119d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1120e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1121d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1122d5b485abSSatish Balay ncols = nzA + nzB; 1123d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1124a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 1125d5b485abSSatish Balay 1126d5b485abSSatish Balay /* load the column values for this row into vals*/ 11275b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 1128d5b485abSSatish Balay for (l=0; l<nzB; l++) { 11293a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 11303eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 11313a40ed3dSBarry Smith } 1132d5b485abSSatish Balay else break; 1133d5b485abSSatish Balay } 1134d5b485abSSatish Balay imark = l; 11353a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 11363eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 11373a40ed3dSBarry Smith } 11383a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 11393eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 11403a40ed3dSBarry Smith } 1141d5b485abSSatish Balay ct2 += ncols; 1142d5b485abSSatish Balay } 1143d5b485abSSatish Balay } 11443eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 1145d5b485abSSatish Balay } 1146d5b485abSSatish Balay } 1147b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 1148b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 1149606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1150d5b485abSSatish Balay 1151d5b485abSSatish Balay /* Form the matrix */ 1152d5b485abSSatish Balay /* create col map */ 1153d5b485abSSatish Balay { 1154d5b485abSSatish Balay int *icol_i; 1155233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1156233c2ff6SSatish Balay /* Create row map*/ 1157b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr); 1158233c2ff6SSatish Balay for (i=0; i<ismax+1; i++) { 1159233c2ff6SSatish Balay ierr = PetscTableCreate(((i<ismax) ? ncol[i] : ncol[i-1])+1,&colmaps[i]);CHKERRQ(ierr); 1160233c2ff6SSatish Balay } 1161233c2ff6SSatish Balay #else 1162a2feddc5SSatish Balay len = (1+ismax)*sizeof(int*)+ ismax*c->Nbs*sizeof(int); 116382502324SSatish Balay ierr = PetscMalloc(len,&cmap);CHKERRQ(ierr); 1164d5b485abSSatish Balay cmap[0] = (int *)(cmap + ismax); 1165549d3d68SSatish Balay ierr = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int));CHKERRQ(ierr); 1166a2feddc5SSatish Balay for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; } 1167233c2ff6SSatish Balay #endif 1168d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1169d5b485abSSatish Balay jmax = ncol[i]; 1170d5b485abSSatish Balay icol_i = icol[i]; 1171233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1172233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1173233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1174001aedefSBarry Smith ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr); 1175233c2ff6SSatish Balay } 1176233c2ff6SSatish Balay #else 1177d5b485abSSatish Balay cmap_i = cmap[i]; 1178d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1179d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1180d5b485abSSatish Balay } 1181233c2ff6SSatish Balay #endif 1182d5b485abSSatish Balay } 1183d5b485abSSatish Balay } 1184d5b485abSSatish Balay 1185d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1186d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1187d5b485abSSatish Balay len = (1+ismax)*sizeof(int*)+ j*sizeof(int); 118882502324SSatish Balay ierr = PetscMalloc(len,&lens);CHKERRQ(ierr); 1189d5b485abSSatish Balay lens[0] = (int *)(lens + ismax); 1190549d3d68SSatish Balay ierr = PetscMemzero(lens[0],j*sizeof(int));CHKERRQ(ierr); 1191d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1192d5b485abSSatish Balay 1193d5b485abSSatish Balay /* Update lens from local data */ 1194d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1195d5b485abSSatish Balay jmax = nrow[i]; 1196233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1197233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1198233c2ff6SSatish Balay #else 1199d5b485abSSatish Balay cmap_i = cmap[i]; 1200233c2ff6SSatish Balay #endif 1201d5b485abSSatish Balay irow_i = irow[i]; 1202d5b485abSSatish Balay lens_i = lens[i]; 1203d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1204d5b485abSSatish Balay row = irow_i[j]; 1205233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1206233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1207233c2ff6SSatish Balay #else 1208d5b485abSSatish Balay proc = rtable[row]; 1209233c2ff6SSatish Balay #endif 1210d5b485abSSatish Balay if (proc == rank) { 121136f4e84dSSatish Balay /* Get indices from matA and then from matB */ 121236f4e84dSSatish Balay row = row - rstart; 121336f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 121436f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1215233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1216233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 1217233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1218233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1219233c2ff6SSatish Balay } 1220233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 1221233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1222233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1223233c2ff6SSatish Balay } 1224233c2ff6SSatish Balay #else 1225ca161407SBarry Smith for (k=0; k<nzA; k++) { 122636f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; } 1227ca161407SBarry Smith } 1228ca161407SBarry Smith for (k=0; k<nzB; k++) { 122936f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; } 1230d5b485abSSatish Balay } 1231233c2ff6SSatish Balay #endif 1232d5b485abSSatish Balay } 1233a2feddc5SSatish Balay } 1234ca161407SBarry Smith } 1235233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1236d5b485abSSatish Balay /* Create row map*/ 1237b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr); 1238233c2ff6SSatish Balay for (i=0; i<ismax+1; i++){ 1239233c2ff6SSatish Balay ierr = PetscTableCreate((i<ismax) ? nrow[i] : nrow[i-1],&rowmaps[i]);CHKERRQ(ierr); 1240233c2ff6SSatish Balay } 1241233c2ff6SSatish Balay #else 1242233c2ff6SSatish Balay /* Create row map*/ 1243233c2ff6SSatish Balay len = (1+ismax)*sizeof(int*)+ ismax*Mbs*sizeof(int); 124482502324SSatish Balay ierr = PetscMalloc(len,&rmap);CHKERRQ(ierr); 1245d5b485abSSatish Balay rmap[0] = (int *)(rmap + ismax); 1246233c2ff6SSatish Balay ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(int));CHKERRQ(ierr); 1247233c2ff6SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;} 1248233c2ff6SSatish Balay #endif 1249d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1250d5b485abSSatish Balay irow_i = irow[i]; 1251d5b485abSSatish Balay jmax = nrow[i]; 1252233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1253233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1254233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1255233c2ff6SSatish Balay ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr); 1256233c2ff6SSatish Balay } 1257233c2ff6SSatish Balay #else 1258233c2ff6SSatish Balay rmap_i = rmap[i]; 1259d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1260d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1261d5b485abSSatish Balay } 1262233c2ff6SSatish Balay #endif 1263d5b485abSSatish Balay } 1264d5b485abSSatish Balay 1265d5b485abSSatish Balay /* Update lens from offproc data */ 1266d5b485abSSatish Balay { 1267d5b485abSSatish Balay int *rbuf2_i,*rbuf3_i,*sbuf1_i; 1268d5b485abSSatish Balay 1269d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1270ca161407SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);CHKERRQ(ierr); 1271d5b485abSSatish Balay index = pa[i]; 1272d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1273d5b485abSSatish Balay jmax = sbuf1_i[0]; 1274d5b485abSSatish Balay ct1 = 2*jmax+1; 1275d5b485abSSatish Balay ct2 = 0; 1276d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1277d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1278d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1279d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1280d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1281d5b485abSSatish Balay lens_i = lens[is_no]; 1282233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1283233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1284233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1285233c2ff6SSatish Balay #else 1286d5b485abSSatish Balay cmap_i = cmap[is_no]; 1287d5b485abSSatish Balay rmap_i = rmap[is_no]; 1288233c2ff6SSatish Balay #endif 1289d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1290233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1291233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1292233c2ff6SSatish Balay row--; 129329bbc08cSBarry Smith if(row < 0) { SETERRQ(1,"row not found in table"); } 1294233c2ff6SSatish Balay #else 1295d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1296233c2ff6SSatish Balay #endif 1297d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1298d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1299233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1300233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 1301233c2ff6SSatish Balay if (tt) { 1302233c2ff6SSatish Balay lens_i[row]++; 1303233c2ff6SSatish Balay } 1304233c2ff6SSatish Balay #else 1305d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1306d5b485abSSatish Balay lens_i[row]++; 1307d5b485abSSatish Balay } 1308233c2ff6SSatish Balay #endif 1309d5b485abSSatish Balay } 1310d5b485abSSatish Balay } 1311d5b485abSSatish Balay } 1312d5b485abSSatish Balay } 1313d5b485abSSatish Balay } 1314606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1315606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1316ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr); 1317606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1318606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1319d5b485abSSatish Balay 1320d5b485abSSatish Balay /* Create the submatrices */ 1321d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1322d5b485abSSatish Balay /* 1323d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1324d5b485abSSatish Balay */ 1325d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1326df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 132736f4e84dSSatish Balay if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) { 132829bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1329d5b485abSSatish Balay } 13300f5bd95cSBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(int),&flag);CHKERRQ(ierr); 13310f5bd95cSBarry Smith if (flag == PETSC_FALSE) { 133229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1333d5b485abSSatish Balay } 1334d5b485abSSatish Balay /* Initial matrix as if empty */ 1335549d3d68SSatish Balay ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(int));CHKERRQ(ierr); 1336ce60de66SLois Curfman McInnes submats[i]->factor = C->factor; 1337d5b485abSSatish Balay } 1338ca161407SBarry Smith } else { 1339d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1340ca161407SBarry Smith ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,a->bs,nrow[i]*bs,ncol[i]*bs,0,lens[i],submats+i);CHKERRQ(ierr); 1341d5b485abSSatish Balay } 1342d5b485abSSatish Balay } 1343d5b485abSSatish Balay 1344d5b485abSSatish Balay /* Assemble the matrices */ 1345d5b485abSSatish Balay /* First assemble the local rows */ 1346d5b485abSSatish Balay { 134736f4e84dSSatish Balay int ilen_row,*imat_ilen,*imat_j,*imat_i; 13483eda8832SBarry Smith MatScalar *imat_a; 1349d5b485abSSatish Balay 1350d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1351df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1352d5b485abSSatish Balay imat_ilen = mat->ilen; 1353d5b485abSSatish Balay imat_j = mat->j; 1354d5b485abSSatish Balay imat_i = mat->i; 1355d5b485abSSatish Balay imat_a = mat->a; 1356233c2ff6SSatish Balay 1357233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1358233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1359233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1360233c2ff6SSatish Balay #else 1361d5b485abSSatish Balay cmap_i = cmap[i]; 1362d5b485abSSatish Balay rmap_i = rmap[i]; 1363233c2ff6SSatish Balay #endif 1364d5b485abSSatish Balay irow_i = irow[i]; 1365d5b485abSSatish Balay jmax = nrow[i]; 1366d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1367d5b485abSSatish Balay row = irow_i[j]; 1368233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1369233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1370233c2ff6SSatish Balay #else 1371d5b485abSSatish Balay proc = rtable[row]; 1372233c2ff6SSatish Balay #endif 1373d5b485abSSatish Balay if (proc == rank) { 137436f4e84dSSatish Balay row = row - rstart; 137536f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 137636f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 137736f4e84dSSatish Balay cworkA = a_j + a_i[row]; 137836f4e84dSSatish Balay cworkB = b_j + b_i[row]; 137936f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 138036f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 1381233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1382233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr); 1383233c2ff6SSatish Balay row--; 138429bbc08cSBarry Smith if (row < 0) { SETERRQ(1,"row not found in table"); } 1385233c2ff6SSatish Balay #else 138636f4e84dSSatish Balay row = rmap_i[row + rstart]; 1387233c2ff6SSatish Balay #endif 1388df36731bSBarry Smith mat_i = imat_i[row]; 138936f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1390d5b485abSSatish Balay mat_j = imat_j + mat_i; 139136f4e84dSSatish Balay ilen_row = imat_ilen[row]; 139236f4e84dSSatish Balay 139336f4e84dSSatish Balay /* load the column indices for this row into cols*/ 139436f4e84dSSatish Balay for (l=0; l<nzB; l++) { 139536f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1396233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1397233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr); 1398233c2ff6SSatish Balay if (tcol) { 1399233c2ff6SSatish Balay #else 140036f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1401233c2ff6SSatish Balay #endif 1402df36731bSBarry Smith *mat_j++ = tcol - 1; 14033eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1404549d3d68SSatish Balay mat_a += bs2; 1405d5b485abSSatish Balay ilen_row++; 1406d5b485abSSatish Balay } 1407ca161407SBarry Smith } else break; 140836f4e84dSSatish Balay } 140936f4e84dSSatish Balay imark = l; 141036f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1411233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1412233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1413233c2ff6SSatish Balay if (tcol) { 1414233c2ff6SSatish Balay #else 141536f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1416233c2ff6SSatish Balay #endif 141736f4e84dSSatish Balay *mat_j++ = tcol - 1; 14183eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1419549d3d68SSatish Balay mat_a += bs2; 142036f4e84dSSatish Balay ilen_row++; 142136f4e84dSSatish Balay } 142236f4e84dSSatish Balay } 142336f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1424233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1425233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1426233c2ff6SSatish Balay if (tcol) { 1427233c2ff6SSatish Balay #else 142836f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1429233c2ff6SSatish Balay #endif 143036f4e84dSSatish Balay *mat_j++ = tcol - 1; 14313eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1432549d3d68SSatish Balay mat_a += bs2; 143336f4e84dSSatish Balay ilen_row++; 143436f4e84dSSatish Balay } 143536f4e84dSSatish Balay } 1436d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1437d5b485abSSatish Balay } 1438d5b485abSSatish Balay } 143936f4e84dSSatish Balay 1440d5b485abSSatish Balay } 1441d5b485abSSatish Balay } 1442d5b485abSSatish Balay 1443d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1444d5b485abSSatish Balay { 1445d5b485abSSatish Balay int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1446d5b485abSSatish Balay int *imat_j,*imat_i; 14473eda8832SBarry Smith MatScalar *imat_a,*rbuf4_i; 1448d5b485abSSatish Balay 1449d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1450ca161407SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);CHKERRQ(ierr); 1451d5b485abSSatish Balay index = pa[i]; 1452d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1453d5b485abSSatish Balay jmax = sbuf1_i[0]; 1454d5b485abSSatish Balay ct1 = 2*jmax + 1; 1455d5b485abSSatish Balay ct2 = 0; 1456d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1457d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1458d5b485abSSatish Balay rbuf4_i = rbuf4[i]; 1459d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1460d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1461233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1462233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1463233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1464233c2ff6SSatish Balay #else 1465d5b485abSSatish Balay rmap_i = rmap[is_no]; 1466d5b485abSSatish Balay cmap_i = cmap[is_no]; 1467233c2ff6SSatish Balay #endif 1468df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1469d5b485abSSatish Balay imat_ilen = mat->ilen; 1470d5b485abSSatish Balay imat_j = mat->j; 1471d5b485abSSatish Balay imat_i = mat->i; 1472d5b485abSSatish Balay imat_a = mat->a; 1473d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1474d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1475d5b485abSSatish Balay row = sbuf1_i[ct1]; 1476233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1477233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr); 1478233c2ff6SSatish Balay row--; 147929bbc08cSBarry Smith if(row < 0) { SETERRQ(1,"row not found in table"); } 1480233c2ff6SSatish Balay #else 1481d5b485abSSatish Balay row = rmap_i[row]; 1482233c2ff6SSatish Balay #endif 1483d5b485abSSatish Balay ilen = imat_ilen[row]; 1484df36731bSBarry Smith mat_i = imat_i[row]; 148536f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1486d5b485abSSatish Balay mat_j = imat_j + mat_i; 1487d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1488d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1489233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1490233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1491233c2ff6SSatish Balay if (tcol) { 1492233c2ff6SSatish Balay #else 1493d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1494233c2ff6SSatish Balay #endif 1495df36731bSBarry Smith *mat_j++ = tcol - 1; 149636f4e84dSSatish Balay /* *mat_a++= rbuf4_i[ct2]; */ 14973eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1498549d3d68SSatish Balay mat_a += bs2; 1499d5b485abSSatish Balay ilen++; 1500d5b485abSSatish Balay } 1501d5b485abSSatish Balay } 1502d5b485abSSatish Balay imat_ilen[row] = ilen; 1503d5b485abSSatish Balay } 1504d5b485abSSatish Balay } 1505d5b485abSSatish Balay } 1506d5b485abSSatish Balay } 1507606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1508606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1509ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr); 1510606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1511606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 1512d5b485abSSatish Balay 1513d5b485abSSatish Balay /* Restore the indices */ 1514d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1515d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1516d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1517d5b485abSSatish Balay } 1518d5b485abSSatish Balay 1519d5b485abSSatish Balay /* Destroy allocated memory */ 1520606d414cSSatish Balay ierr = PetscFree(irow);CHKERRQ(ierr); 1521606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 1522606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1523d5b485abSSatish Balay 1524606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1525606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1526d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1527606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1528d5b485abSSatish Balay } 1529d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1530606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1531606d414cSSatish Balay ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1532d5b485abSSatish Balay } 1533d5b485abSSatish Balay 1534606d414cSSatish Balay ierr = PetscFree(sbuf2);CHKERRQ(ierr); 1535606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1536606d414cSSatish Balay ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1537606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1538606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1539606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1540606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1541d5b485abSSatish Balay 1542233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1543233c2ff6SSatish Balay for (i=0; i<ismax+1; i++){ 1544233c2ff6SSatish Balay ierr = PetscTableDelete(rowmaps[i]);CHKERRQ(ierr); 1545233c2ff6SSatish Balay ierr = PetscTableDelete(colmaps[i]);CHKERRQ(ierr); 1546233c2ff6SSatish Balay } 1547233c2ff6SSatish Balay ierr = PetscFree(colmaps);CHKERRQ(ierr); 1548233c2ff6SSatish Balay ierr = PetscFree(rowmaps);CHKERRQ(ierr); 1549233c2ff6SSatish Balay /* Mark Adams */ 1550233c2ff6SSatish Balay #else 1551606d414cSSatish Balay ierr = PetscFree(rmap);CHKERRQ(ierr); 1552233c2ff6SSatish Balay ierr = PetscFree(cmap);CHKERRQ(ierr); 1553233c2ff6SSatish Balay #endif 1554606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1555d5b485abSSatish Balay 1556d5b485abSSatish Balay for (i=0; i<ismax; i++) { 155736f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 155836f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1559d5b485abSSatish Balay } 156034e6ae18SSatish Balay 15613a40ed3dSBarry Smith PetscFunctionReturn(0); 1562d5b485abSSatish Balay } 1563ca161407SBarry Smith 1564