173f4d377SMatthew Knepley /*$Id: baijov.c,v 1.65 2001/08/06 21:15:42 bsmith Exp $*/ 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*); 1387828ca2SBarry Smith EXTERN int MatGetRow_MPIBAIJ(Mat,int,int*,int**,PetscScalar**); 1487828ca2SBarry Smith EXTERN int MatRestoreRow_MPIBAIJ(Mat,int,int*,int**,PetscScalar**); 15d5b485abSSatish Balay 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "MatCompressIndicesGeneral_MPIBAIJ" 18268466fbSBarry Smith static int MatCompressIndicesGeneral_MPIBAIJ(Mat C,int imax,const 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 884a2ae208SSatish Balay #undef __FUNCT__ 894a2ae208SSatish Balay #define __FUNCT__ "MatCompressIndicesSorted_MPIBAIJ" 90268466fbSBarry Smith static int MatCompressIndicesSorted_MPIBAIJ(Mat C,int imax,const 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 1434a2ae208SSatish Balay #undef __FUNCT__ 1444a2ae208SSatish Balay #define __FUNCT__ "MatExpandIndices_MPIBAIJ" 145268466fbSBarry Smith static int MatExpandIndices_MPIBAIJ(Mat C,int imax,const 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 154268466fbSBarry Smith 1553a40ed3dSBarry Smith PetscFunctionBegin; 156df36731bSBarry Smith 157233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 158233c2ff6SSatish Balay /* Now check max size */ 159233c2ff6SSatish Balay for (i=0,maxsz=0; i<imax; i++) { 160233c2ff6SSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 161233c2ff6SSatish Balay ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 162233c2ff6SSatish Balay if (n*bs > maxsz) maxsz = n*bs; 163233c2ff6SSatish Balay } 164b0a32e0cSBarry Smith ierr = PetscMalloc((maxsz+1)*sizeof(int),&nidx);CHKERRQ(ierr); 165233c2ff6SSatish Balay #else 166b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs*bs+1)*sizeof(int),&nidx);CHKERRQ(ierr); 167233c2ff6SSatish Balay #endif 168df36731bSBarry Smith 169df36731bSBarry Smith for (i=0; i<imax; i++) { 17036f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx);CHKERRQ(ierr); 171b9b97703SBarry Smith ierr = ISGetLocalSize(is_in[i],&n);CHKERRQ(ierr); 172df36731bSBarry Smith for (j=0; j<n ; ++j){ 173df36731bSBarry Smith for (k=0; k<bs; k++) 174df36731bSBarry Smith nidx[j*bs+k] = idx[j]*bs+k; 175df36731bSBarry Smith } 17636f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx);CHKERRQ(ierr); 177329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n*bs,nidx,is_out+i);CHKERRQ(ierr); 178df36731bSBarry Smith } 179606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 1803a40ed3dSBarry Smith PetscFunctionReturn(0); 181df36731bSBarry Smith } 182df36731bSBarry Smith 183df36731bSBarry Smith 1844a2ae208SSatish Balay #undef __FUNCT__ 1854a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ" 186268466fbSBarry Smith int MatIncreaseOverlap_MPIBAIJ(Mat C,int imax,IS is[],int ov) 187d5b485abSSatish Balay { 188d5b485abSSatish Balay int i,ierr; 18936f4e84dSSatish Balay IS *is_new; 19036f4e84dSSatish Balay 1913a40ed3dSBarry Smith PetscFunctionBegin; 19282502324SSatish Balay ierr = PetscMalloc(imax*sizeof(IS),&is_new);CHKERRQ(ierr); 193df36731bSBarry Smith /* Convert the indices into block format */ 194e757655aSSatish Balay ierr = MatCompressIndicesGeneral_MPIBAIJ(C,imax,is,is_new);CHKERRQ(ierr); 19529bbc08cSBarry Smith if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");} 196d5b485abSSatish Balay for (i=0; i<ov; ++i) { 19736f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 198d5b485abSSatish Balay } 199ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);} 20036f4e84dSSatish Balay ierr = MatExpandIndices_MPIBAIJ(C,imax,is_new,is);CHKERRQ(ierr); 201ca161407SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);} 202606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 2033a40ed3dSBarry Smith PetscFunctionReturn(0); 204d5b485abSSatish Balay } 205d5b485abSSatish Balay 206d5b485abSSatish Balay /* 207d5b485abSSatish Balay Sample message format: 208d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 209*0e9b0e7eSHong Zhang to index sets is[1], is[5] 210d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 211d5b485abSSatish Balay ----------- 212d5b485abSSatish Balay mesg [1] = 1 => is[1] 213d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 214d5b485abSSatish Balay ----------- 215d5b485abSSatish Balay mesg [5] = 5 => is[5] 216d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 217d5b485abSSatish Balay ----------- 218d5b485abSSatish Balay mesg [7] 219*0e9b0e7eSHong Zhang mesg [n] data(is[1]) 220d5b485abSSatish Balay ----------- 221d5b485abSSatish Balay mesg[n+1] 222d5b485abSSatish Balay mesg[m] data(is[5]) 223d5b485abSSatish Balay ----------- 224d5b485abSSatish Balay 225d5b485abSSatish Balay Notes: 226d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 227d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 228d5b485abSSatish Balay */ 2294a2ae208SSatish Balay #undef __FUNCT__ 2304a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Once" 231268466fbSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat C,int imax,IS is[]) 232d5b485abSSatish Balay { 233df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 234d5b485abSSatish Balay int **idx,*n,*w1,*w2,*w3,*w4,*rtable,**data,len,*idx_i; 235df36731bSBarry Smith int size,rank,Mbs,i,j,k,ierr,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 236c7dd2924SSatish Balay int *ctr,*pa,*tmp,nrqr,*isz,*isz1,**xdata,**rbuf2; 237c7dd2924SSatish Balay int *onodes1,*olengths1,tag1,tag2,*onodes2,*olengths2; 238f1af5d2fSBarry Smith PetscBT *table; 239d5b485abSSatish Balay MPI_Comm comm; 240d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 241d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 242d5b485abSSatish Balay 2433a40ed3dSBarry Smith PetscFunctionBegin; 244d5b485abSSatish Balay comm = C->comm; 245d5b485abSSatish Balay size = c->size; 246d5b485abSSatish Balay rank = c->rank; 247df36731bSBarry Smith Mbs = c->Mbs; 248d5b485abSSatish Balay 249c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 250c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 251c7dd2924SSatish Balay 252df36731bSBarry Smith len = (imax+1)*sizeof(int*)+ (imax + Mbs)*sizeof(int); 25382502324SSatish Balay ierr = PetscMalloc(len,&idx);CHKERRQ(ierr); 254d5b485abSSatish Balay n = (int*)(idx + imax); 255d5b485abSSatish Balay rtable = n + imax; 256d5b485abSSatish Balay 257d5b485abSSatish Balay for (i=0; i<imax; i++) { 258d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 259b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 260d5b485abSSatish Balay } 261d5b485abSSatish Balay 262d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 263d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 264d5b485abSSatish Balay len = c->rowners[i+1]; 265d5b485abSSatish Balay for (; j<len; j++) { 266d5b485abSSatish Balay rtable[j] = i; 267d5b485abSSatish Balay } 268d5b485abSSatish Balay } 269d5b485abSSatish Balay 270d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 271d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 27282502324SSatish Balay ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr);/* mesg size */ 273d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 274d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 275d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 276549d3d68SSatish Balay ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 277d5b485abSSatish Balay for (i=0; i<imax; i++) { 278549d3d68SSatish Balay ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 279d5b485abSSatish Balay idx_i = idx[i]; 280d5b485abSSatish Balay len = n[i]; 281d5b485abSSatish Balay for (j=0; j<len; j++) { 282d5b485abSSatish Balay row = idx_i[j]; 2836b41c64dSBarry Smith if (row < 0) { 28429bbc08cSBarry Smith SETERRQ(1,"Index set cannot have negative entries"); 2856b41c64dSBarry Smith } 286d5b485abSSatish Balay proc = rtable[row]; 287d5b485abSSatish Balay w4[proc]++; 288d5b485abSSatish Balay } 289d5b485abSSatish Balay for (j=0; j<size; j++){ 290d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 291d5b485abSSatish Balay } 292d5b485abSSatish Balay } 293d5b485abSSatish Balay 294d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 295d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 296*0e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 297d5b485abSSatish Balay w3[rank] = 0; 298d5b485abSSatish Balay for (i=0; i<size; i++) { 299d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 300d5b485abSSatish Balay } 301d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 302b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); 303d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 304d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 305d5b485abSSatish Balay } 306d5b485abSSatish Balay 307d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 308d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 309d5b485abSSatish Balay j = pa[i]; 310d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 311d5b485abSSatish Balay msz += w1[j]; 312d5b485abSSatish Balay } 313d5b485abSSatish Balay 314c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 315a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 316c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 317d5b485abSSatish Balay 318c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 319c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 320d5b485abSSatish Balay 321d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 322d5b485abSSatish Balay len = 2*size*sizeof(int*) + (size+msz)*sizeof(int); 32382502324SSatish Balay ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr); 324d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 325549d3d68SSatish Balay ierr = PetscMemzero(outdat,2*size*sizeof(int*));CHKERRQ(ierr); 326d5b485abSSatish Balay tmp = (int*)(outdat + 2*size); 327d5b485abSSatish Balay ctr = tmp + msz; 328d5b485abSSatish Balay 329d5b485abSSatish Balay { 330d5b485abSSatish Balay int *iptr = tmp,ict = 0; 331d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 332d5b485abSSatish Balay j = pa[i]; 333d5b485abSSatish Balay iptr += ict; 334d5b485abSSatish Balay outdat[j] = iptr; 335d5b485abSSatish Balay ict = w1[j]; 336d5b485abSSatish Balay } 337d5b485abSSatish Balay } 338d5b485abSSatish Balay 339d5b485abSSatish Balay /* Form the outgoing messages */ 340d5b485abSSatish Balay /*plug in the headers*/ 341d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 342d5b485abSSatish Balay j = pa[i]; 343d5b485abSSatish Balay outdat[j][0] = 0; 344549d3d68SSatish Balay ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 345d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 346d5b485abSSatish Balay } 347d5b485abSSatish Balay 348d5b485abSSatish Balay /* Memory for doing local proc's work*/ 349d5b485abSSatish Balay { 350d5b485abSSatish Balay int *d_p; 351d5b485abSSatish Balay char *t_p; 352d5b485abSSatish Balay 353f1af5d2fSBarry Smith len = (imax)*(sizeof(PetscBT) + sizeof(int*)+ sizeof(int)) + 354b6410449SSatish Balay (Mbs)*imax*sizeof(int) + (Mbs/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char) + 1; 35582502324SSatish Balay ierr = PetscMalloc(len,&table);CHKERRQ(ierr); 356549d3d68SSatish Balay ierr = PetscMemzero(table,len);CHKERRQ(ierr); 357d5b485abSSatish Balay data = (int **)(table + imax); 358bbd702dbSSatish Balay isz = (int *)(data + imax); 359bbd702dbSSatish Balay d_p = (int *)(isz + imax); 360bbd702dbSSatish Balay t_p = (char *)(d_p + Mbs*imax); 361bbd702dbSSatish Balay for (i=0; i<imax; i++) { 362b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 363bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 364d5b485abSSatish Balay } 365d5b485abSSatish Balay } 366d5b485abSSatish Balay 367d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 368d5b485abSSatish Balay { 369d5b485abSSatish Balay int n_i,*data_i,isz_i,*outdat_j,ctr_j; 370f1af5d2fSBarry Smith PetscBT table_i; 371d5b485abSSatish Balay 372d5b485abSSatish Balay for (i=0; i<imax; i++) { 373549d3d68SSatish Balay ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 374d5b485abSSatish Balay n_i = n[i]; 375d5b485abSSatish Balay table_i = table[i]; 376d5b485abSSatish Balay idx_i = idx[i]; 377d5b485abSSatish Balay data_i = data[i]; 378d5b485abSSatish Balay isz_i = isz[i]; 379d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 380d5b485abSSatish Balay row = idx_i[j]; 381d5b485abSSatish Balay proc = rtable[row]; 382d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 383d5b485abSSatish Balay ctr[proc]++; 384d5b485abSSatish Balay *ptr[proc] = row; 385d5b485abSSatish Balay ptr[proc]++; 386d5b485abSSatish Balay } 387d5b485abSSatish Balay else { /* Update the local table */ 388f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 389d5b485abSSatish Balay } 390d5b485abSSatish Balay } 391d5b485abSSatish Balay /* Update the headers for the current IS */ 392d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 393d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 394d5b485abSSatish Balay outdat_j = outdat[j]; 395d5b485abSSatish Balay k = ++outdat_j[0]; 396d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 397d5b485abSSatish Balay outdat_j[2*k-1] = i; 398d5b485abSSatish Balay } 399d5b485abSSatish Balay } 400d5b485abSSatish Balay isz[i] = isz_i; 401d5b485abSSatish Balay } 402d5b485abSSatish Balay } 403d5b485abSSatish Balay 404d5b485abSSatish Balay /* Now post the sends */ 405b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 406d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 407d5b485abSSatish Balay j = pa[i]; 408c7dd2924SSatish Balay ierr = MPI_Isend(outdat[j],w1[j],MPI_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 409d5b485abSSatish Balay } 410d5b485abSSatish Balay 411d5b485abSSatish Balay /* No longer need the original indices*/ 412d5b485abSSatish Balay for (i=0; i<imax; ++i) { 413d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 414d5b485abSSatish Balay } 415606d414cSSatish Balay ierr = PetscFree(idx);CHKERRQ(ierr); 416d5b485abSSatish Balay 417d5b485abSSatish Balay for (i=0; i<imax; ++i) { 418d5b485abSSatish Balay ierr = ISDestroy(is[i]);CHKERRQ(ierr); 419d5b485abSSatish Balay } 420d5b485abSSatish Balay 421d5b485abSSatish Balay /* Do Local work*/ 422df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 423d5b485abSSatish Balay 424d5b485abSSatish Balay /* Receive messages*/ 425b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 426c7dd2924SSatish Balay ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr); 427d5b485abSSatish Balay 428b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 429ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr); 430d5b485abSSatish Balay 431d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 432606d414cSSatish Balay ierr = PetscFree(outdat);CHKERRQ(ierr); 433606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 434d5b485abSSatish Balay 43582502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(int *),&xdata);CHKERRQ(ierr); 436b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(int),&isz1);CHKERRQ(ierr); 437df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 438606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 439d5b485abSSatish Balay 440d5b485abSSatish Balay /* Send the data back*/ 441d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 442d5b485abSSatish Balay { 443c7dd2924SSatish Balay int *rw1; 444d5b485abSSatish Balay 445c7dd2924SSatish Balay ierr = PetscMalloc(size*sizeof(int),&rw1);CHKERRQ(ierr); 446c7dd2924SSatish Balay ierr = PetscMemzero(rw1,size*sizeof(int));CHKERRQ(ierr); 447c7dd2924SSatish Balay 448d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 449d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 450c7dd2924SSatish Balay if (proc != onodes1[i]) SETERRQ(1,"MPI_SOURCE mismatch"); 451d5b485abSSatish Balay rw1[proc] = isz1[i]; 452d5b485abSSatish Balay } 453d5b485abSSatish Balay 454c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 455c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 456d5b485abSSatish Balay 457c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 458c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 459c7dd2924SSatish Balay PetscFree(rw1); 460c7dd2924SSatish Balay } 461c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 462c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 463d5b485abSSatish Balay 464d5b485abSSatish Balay /* Now post the sends */ 465b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 466d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 467d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 468c7dd2924SSatish Balay ierr = MPI_Isend(xdata[i],isz1[i],MPI_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 469d5b485abSSatish Balay } 470d5b485abSSatish Balay 471d5b485abSSatish Balay /* receive work done on other processors*/ 472d5b485abSSatish Balay { 473999d9058SBarry Smith int idex,is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 474f1af5d2fSBarry Smith PetscBT table_i; 475d5b485abSSatish Balay MPI_Status *status2; 476d5b485abSSatish Balay 477169449a0SSatish Balay ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 478d5b485abSSatish Balay 479d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 480999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 481d5b485abSSatish Balay /* Process the message*/ 482999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 483d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 484999d9058SBarry Smith jmax = rbuf2[idex][0]; 485d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 486d5b485abSSatish Balay max = rbuf2_i[2*j]; 487d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 488d5b485abSSatish Balay isz_i = isz[is_no]; 489d5b485abSSatish Balay data_i = data[is_no]; 490d5b485abSSatish Balay table_i = table[is_no]; 491d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 492d5b485abSSatish Balay row = rbuf2_i[ct1]; 493f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 494d5b485abSSatish Balay } 495d5b485abSSatish Balay isz[is_no] = isz_i; 496d5b485abSSatish Balay } 497d5b485abSSatish Balay } 498ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr); 499606d414cSSatish Balay ierr = PetscFree(status2);CHKERRQ(ierr); 500d5b485abSSatish Balay } 501d5b485abSSatish Balay 502d5b485abSSatish Balay for (i=0; i<imax; ++i) { 503029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr); 504d5b485abSSatish Balay } 505d5b485abSSatish Balay 506c7dd2924SSatish Balay 507c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 508c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 509c7dd2924SSatish Balay 510606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 511606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 512606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 513606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 514606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 515606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 516606d414cSSatish Balay ierr = PetscFree(table);CHKERRQ(ierr); 517606d414cSSatish Balay ierr = PetscFree(s_status);CHKERRQ(ierr); 518606d414cSSatish Balay ierr = PetscFree(recv_status);CHKERRQ(ierr); 519606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 520606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 521606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 5223a40ed3dSBarry Smith PetscFunctionReturn(0); 523d5b485abSSatish Balay } 524d5b485abSSatish Balay 5254a2ae208SSatish Balay #undef __FUNCT__ 5264a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Local" 527d5b485abSSatish Balay /* 528df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 529d5b485abSSatish Balay the work on the local processor. 530d5b485abSSatish Balay 531d5b485abSSatish Balay Inputs: 532df36731bSBarry Smith C - MAT_MPIBAIJ; 533d5b485abSSatish Balay imax - total no of index sets processed at a time; 534df36731bSBarry Smith table - an array of char - size = Mbs bits. 535d5b485abSSatish Balay 536d5b485abSSatish Balay Output: 537*0e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 538d5b485abSSatish Balay to each index set; 539d5b485abSSatish Balay data - pointer to the solutions 540d5b485abSSatish Balay */ 541f1af5d2fSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,PetscBT *table,int *isz,int **data) 542d5b485abSSatish Balay { 543df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 544d5b485abSSatish Balay Mat A = c->A,B = c->B; 545df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 546df36731bSBarry Smith int start,end,val,max,rstart,cstart,*ai,*aj; 547d5b485abSSatish Balay int *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 548f1af5d2fSBarry Smith PetscBT table_i; 549d5b485abSSatish Balay 5503a40ed3dSBarry Smith PetscFunctionBegin; 551d5b485abSSatish Balay rstart = c->rstart; 552d5b485abSSatish Balay cstart = c->cstart; 553d5b485abSSatish Balay ai = a->i; 554df36731bSBarry Smith aj = a->j; 555d5b485abSSatish Balay bi = b->i; 556df36731bSBarry Smith bj = b->j; 557d5b485abSSatish Balay garray = c->garray; 558d5b485abSSatish Balay 559d5b485abSSatish Balay 560d5b485abSSatish Balay for (i=0; i<imax; i++) { 561d5b485abSSatish Balay data_i = data[i]; 562d5b485abSSatish Balay table_i = table[i]; 563d5b485abSSatish Balay isz_i = isz[i]; 564d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 565d5b485abSSatish Balay row = data_i[j] - rstart; 566d5b485abSSatish Balay start = ai[row]; 567d5b485abSSatish Balay end = ai[row+1]; 568d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 569df36731bSBarry Smith val = aj[k] + cstart; 570f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 571d5b485abSSatish Balay } 572d5b485abSSatish Balay start = bi[row]; 573d5b485abSSatish Balay end = bi[row+1]; 574d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 575df36731bSBarry Smith val = garray[bj[k]]; 576f1af5d2fSBarry Smith if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 577d5b485abSSatish Balay } 578d5b485abSSatish Balay } 579d5b485abSSatish Balay isz[i] = isz_i; 580d5b485abSSatish Balay } 5813a40ed3dSBarry Smith PetscFunctionReturn(0); 582d5b485abSSatish Balay } 5834a2ae208SSatish Balay #undef __FUNCT__ 5844a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_MPIBAIJ_Receive" 585d5b485abSSatish Balay /* 586df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 587d5b485abSSatish Balay and return the output 588d5b485abSSatish Balay 589d5b485abSSatish Balay Input: 590d5b485abSSatish Balay C - the matrix 591d5b485abSSatish Balay nrqr - no of messages being processed. 592d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 593d5b485abSSatish Balay 594d5b485abSSatish Balay Output: 595d5b485abSSatish Balay xdata - array of messages to be sent back 596d5b485abSSatish Balay isz1 - size of each message 597d5b485abSSatish Balay 598d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 599d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 600*0e9b0e7eSHong Zhang rather than all previous rows as it is now where a single large chunck of 601d5b485abSSatish Balay memory is used. 602d5b485abSSatish Balay 603d5b485abSSatish Balay */ 604701b8814SBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf,int **xdata,int * isz1) 605d5b485abSSatish Balay { 606df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 607d5b485abSSatish Balay Mat A = c->A,B = c->B; 608df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 609df36731bSBarry Smith int rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 610d5b485abSSatish Balay int row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 611df36731bSBarry Smith int val,max1,max2,rank,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 612bbd702dbSSatish Balay int *rbuf_i,kmax,rbuf_0,ierr; 613f1af5d2fSBarry Smith PetscBT xtable; 614d5b485abSSatish Balay 6153a40ed3dSBarry Smith PetscFunctionBegin; 616d5b485abSSatish Balay rank = c->rank; 617df36731bSBarry Smith Mbs = c->Mbs; 618d5b485abSSatish Balay rstart = c->rstart; 619d5b485abSSatish Balay cstart = c->cstart; 620d5b485abSSatish Balay ai = a->i; 621df36731bSBarry Smith aj = a->j; 622d5b485abSSatish Balay bi = b->i; 623df36731bSBarry Smith bj = b->j; 624d5b485abSSatish Balay garray = c->garray; 625d5b485abSSatish Balay 626d5b485abSSatish Balay 627d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 628d5b485abSSatish Balay rbuf_i = rbuf[i]; 629d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 630d5b485abSSatish Balay ct += rbuf_0; 631d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 632d5b485abSSatish Balay } 633d5b485abSSatish Balay 634701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 635701b8814SBarry Smith else max1 = 1; 636d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 63782502324SSatish Balay ierr = PetscMalloc(mem_estimate*sizeof(int),&xdata[0]);CHKERRQ(ierr); 638d5b485abSSatish Balay ++no_malloc; 6396831982aSBarry Smith ierr = PetscBTCreate(Mbs,xtable);CHKERRQ(ierr); 640549d3d68SSatish Balay ierr = PetscMemzero(isz1,nrqr*sizeof(int));CHKERRQ(ierr); 641d5b485abSSatish Balay 642d5b485abSSatish Balay ct3 = 0; 643d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 644d5b485abSSatish Balay rbuf_i = rbuf[i]; 645d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 646d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 647d5b485abSSatish Balay ct2 = ct1; 648d5b485abSSatish Balay ct3 += ct1; 649d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 6506831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 651d5b485abSSatish Balay oct2 = ct2; 652d5b485abSSatish Balay kmax = rbuf_i[2*j]; 653d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 654d5b485abSSatish Balay row = rbuf_i[ct1]; 655f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 656d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 657d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 65882502324SSatish Balay ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr); 659549d3d68SSatish Balay ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 660606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 661d5b485abSSatish Balay xdata[0] = tmp; 662d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 663d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 664d5b485abSSatish Balay } 665d5b485abSSatish Balay xdata[i][ct2++] = row; 666d5b485abSSatish Balay ct3++; 667d5b485abSSatish Balay } 668d5b485abSSatish Balay } 669d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 670d5b485abSSatish Balay row = xdata[i][k] - rstart; 671d5b485abSSatish Balay start = ai[row]; 672d5b485abSSatish Balay end = ai[row+1]; 673d5b485abSSatish Balay for (l=start; l<end; l++) { 674df36731bSBarry Smith val = aj[l] + cstart; 675f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 676d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 677d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 67882502324SSatish Balay ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr); 679549d3d68SSatish Balay ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 680606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 681d5b485abSSatish Balay xdata[0] = tmp; 682d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 683d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 684d5b485abSSatish Balay } 685d5b485abSSatish Balay xdata[i][ct2++] = val; 686d5b485abSSatish Balay ct3++; 687d5b485abSSatish Balay } 688d5b485abSSatish Balay } 689d5b485abSSatish Balay start = bi[row]; 690d5b485abSSatish Balay end = bi[row+1]; 691d5b485abSSatish Balay for (l=start; l<end; l++) { 692df36731bSBarry Smith val = garray[bj[l]]; 693f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 694d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 695d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 69682502324SSatish Balay ierr = PetscMalloc(new_estimate * sizeof(int),&tmp);CHKERRQ(ierr); 697549d3d68SSatish Balay ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 698606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 699d5b485abSSatish Balay xdata[0] = tmp; 700d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 701d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 702d5b485abSSatish Balay } 703d5b485abSSatish Balay xdata[i][ct2++] = val; 704d5b485abSSatish Balay ct3++; 705d5b485abSSatish Balay } 706d5b485abSSatish Balay } 707d5b485abSSatish Balay } 708d5b485abSSatish Balay /* Update the header*/ 709d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 710d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 711d5b485abSSatish Balay } 712d5b485abSSatish Balay xdata[i][0] = rbuf_0; 713d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 714d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 715d5b485abSSatish Balay } 7166831982aSBarry Smith ierr = PetscBTDestroy(xtable);CHKERRQ(ierr); 717b0a32e0cSBarry Smith PetscLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d, no of mallocs = %d\n",rank,mem_estimate,ct3,no_malloc); 7183a40ed3dSBarry Smith PetscFunctionReturn(0); 719d5b485abSSatish Balay } 720d5b485abSSatish Balay 721268466fbSBarry Smith static int MatGetSubMatrices_MPIBAIJ_local(Mat,int,const IS[],const IS[],MatReuse,Mat *); 722a2feddc5SSatish Balay 7234a2ae208SSatish Balay #undef __FUNCT__ 7244a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ" 725268466fbSBarry Smith int MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 726d5b485abSSatish Balay { 72736f4e84dSSatish Balay IS *isrow_new,*iscol_new; 728cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 729cf4f527aSSatish Balay int nmax,nstages_local,nstages,i,pos,max_no,ierr; 730a2feddc5SSatish Balay 7313a40ed3dSBarry Smith PetscFunctionBegin; 7323a40ed3dSBarry Smith /* The compression and expansion should be avoided. Does'nt point 733a2feddc5SSatish Balay out errors might change the indices hence buggey */ 734a2feddc5SSatish Balay 735b0a32e0cSBarry Smith ierr = PetscMalloc(2*(ismax+1)*sizeof(IS),&isrow_new);CHKERRQ(ierr); 73636f4e84dSSatish Balay iscol_new = isrow_new + ismax; 737e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C,ismax,isrow,isrow_new);CHKERRQ(ierr); 738e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C,ismax,iscol,iscol_new);CHKERRQ(ierr); 739cf4f527aSSatish Balay 740cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 741cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 74282502324SSatish Balay ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 743cf4f527aSSatish Balay } 744cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 745cf4f527aSSatish Balay nmax = 20*1000000 / (c->Nbs * sizeof(int)); 746329f5518SBarry Smith if (!nmax) nmax = 1; 747cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 748cf4f527aSSatish Balay 749653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 750ca161407SBarry Smith ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 751cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 752cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 753cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 754cf4f527aSSatish Balay else max_no = ismax-pos; 755cf4f527aSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos);CHKERRQ(ierr); 756cf4f527aSSatish Balay pos += max_no; 757cf4f527aSSatish Balay } 75836f4e84dSSatish Balay 75936f4e84dSSatish Balay for (i=0; i<ismax; i++) { 760ca161407SBarry Smith ierr = ISDestroy(isrow_new[i]);CHKERRQ(ierr); 761ca161407SBarry Smith ierr = ISDestroy(iscol_new[i]);CHKERRQ(ierr); 76236f4e84dSSatish Balay } 763606d414cSSatish Balay ierr = PetscFree(isrow_new);CHKERRQ(ierr); 7643a40ed3dSBarry Smith PetscFunctionReturn(0); 765a2feddc5SSatish Balay } 766a2feddc5SSatish Balay 767233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 7684a2ae208SSatish Balay #undef __FUNCT__ 7694a2ae208SSatish Balay #define __FUNCT__ "PetscGetProc" 770233c2ff6SSatish Balay int PetscGetProc(const int gid, const int numprocs, const int proc_gnode[], int *proc) 771233c2ff6SSatish Balay { 772233c2ff6SSatish Balay int nGlobalNd = proc_gnode[numprocs]; 773233c2ff6SSatish Balay int fproc = (int) ((float)gid * (float)numprocs / (float)nGlobalNd + 0.5); 774233c2ff6SSatish Balay 775233c2ff6SSatish Balay PetscFunctionBegin; 77629bbc08cSBarry Smith /* if(fproc < 0) SETERRQ(1,"fproc < 0");*/ 777233c2ff6SSatish Balay if (fproc > numprocs) fproc = numprocs; 778233c2ff6SSatish Balay while (gid < proc_gnode[fproc] || gid >= proc_gnode[fproc+1]) { 779233c2ff6SSatish Balay if (gid < proc_gnode[fproc]) fproc--; 780233c2ff6SSatish Balay else fproc++; 781233c2ff6SSatish Balay } 78229bbc08cSBarry Smith /* if(fproc<0 || fproc>=numprocs) { SETERRQ(1,"fproc < 0 || fproc >= numprocs"); }*/ 783233c2ff6SSatish Balay *proc = fproc; 784233c2ff6SSatish Balay PetscFunctionReturn(0); 785233c2ff6SSatish Balay } 786233c2ff6SSatish Balay #endif 787233c2ff6SSatish Balay 788a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 7894a2ae208SSatish Balay #undef __FUNCT__ 7904a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_MPIBAIJ_local" 791268466fbSBarry Smith static int MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 792a2feddc5SSatish Balay { 793df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 794cf4f527aSSatish Balay Mat A = c->A; 795df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 796233c2ff6SSatish Balay int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,start,end,size; 797233c2ff6SSatish Balay int **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,ierr,**rbuf1,row,proc; 798999d9058SBarry Smith int nrqs,msz,**ptr,idex,*req_size,*ctr,*pa,*tmp,tcol,nrqr; 799233c2ff6SSatish Balay int **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 800233c2ff6SSatish Balay int **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 801233c2ff6SSatish Balay int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 802233c2ff6SSatish Balay int bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 80336f4e84dSSatish Balay int cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 80434e6ae18SSatish Balay int *bmap = c->garray,ctmp,rstart=c->rstart,tag0,tag1,tag2,tag3; 805d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 806d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 807d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 808d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 809d5b485abSSatish Balay MPI_Comm comm; 8103eda8832SBarry Smith MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 8113eda8832SBarry Smith MatScalar *a_a=a->a,*b_a=b->a; 8120f5bd95cSBarry Smith PetscTruth flag; 813c7dd2924SSatish Balay int *onodes1,*olengths1; 814c7dd2924SSatish Balay 81580d608c0SSatish Balay #if defined (PETSC_USE_CTABLE) 816233c2ff6SSatish Balay int tt; 817233c2ff6SSatish Balay PetscTable *rowmaps,*colmaps,lrow1_grow1,lcol1_gcol1; 818233c2ff6SSatish Balay #else 819233c2ff6SSatish Balay int **cmap,*cmap_i,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 820233c2ff6SSatish Balay #endif 821d5b485abSSatish Balay 8223a40ed3dSBarry Smith PetscFunctionBegin; 823d5b485abSSatish Balay comm = C->comm; 82434e6ae18SSatish Balay tag0 = C->tag; 825d5b485abSSatish Balay size = c->size; 826d5b485abSSatish Balay rank = c->rank; 827d5b485abSSatish Balay 82834e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 82934e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 83034e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 83134e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 83234e6ae18SSatish Balay 833d5b485abSSatish Balay /* Check if the col indices are sorted */ 834d5b485abSSatish Balay for (i=0; i<ismax; i++) { 835888f2ed8SSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr); 83629bbc08cSBarry Smith if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 837d5b485abSSatish Balay } 838d5b485abSSatish Balay 839233c2ff6SSatish Balay len = (2*ismax+1)*(sizeof(int*)+ sizeof(int)); 840233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 841233c2ff6SSatish Balay len += (Mbs+1)*sizeof(int); 842233c2ff6SSatish Balay #endif 84382502324SSatish Balay ierr = PetscMalloc(len,&irow);CHKERRQ(ierr); 844d5b485abSSatish Balay icol = irow + ismax; 845d5b485abSSatish Balay nrow = (int*)(icol + ismax); 846d5b485abSSatish Balay ncol = nrow + ismax; 847233c2ff6SSatish Balay #if !defined (PETSC_USE_CTABLE) 848d5b485abSSatish Balay rtable = ncol + ismax; 849d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 850d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 851d5b485abSSatish Balay jmax = c->rowners[i+1]; 852d5b485abSSatish Balay for (; j<jmax; j++) { 853d5b485abSSatish Balay rtable[j] = i; 854d5b485abSSatish Balay } 855d5b485abSSatish Balay } 856233c2ff6SSatish Balay #endif 857233c2ff6SSatish Balay 858233c2ff6SSatish Balay for (i=0; i<ismax; i++) { 859233c2ff6SSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 860233c2ff6SSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 861233c2ff6SSatish Balay ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 862233c2ff6SSatish Balay ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 863233c2ff6SSatish Balay } 864d5b485abSSatish Balay 865d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg,and buffer space 866d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 86782502324SSatish Balay ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr); /* mesg size */ 868d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 869d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 870d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 871549d3d68SSatish Balay ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 872d5b485abSSatish Balay for (i=0; i<ismax; i++) { 873549d3d68SSatish Balay ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 874d5b485abSSatish Balay jmax = nrow[i]; 875d5b485abSSatish Balay irow_i = irow[i]; 876d5b485abSSatish Balay for (j=0; j<jmax; j++) { 877d5b485abSSatish Balay row = irow_i[j]; 878233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 879233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 880233c2ff6SSatish Balay #else 881d5b485abSSatish Balay proc = rtable[row]; 882233c2ff6SSatish Balay #endif 883d5b485abSSatish Balay w4[proc]++; 884d5b485abSSatish Balay } 885d5b485abSSatish Balay for (j=0; j<size; j++) { 886d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 887d5b485abSSatish Balay } 888d5b485abSSatish Balay } 889d5b485abSSatish Balay 890d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 891e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 892d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 893d5b485abSSatish Balay w3[rank] = 0; 894d5b485abSSatish Balay for (i=0; i<size; i++) { 895d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 896d5b485abSSatish Balay } 897b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); /*(proc -array)*/ 898d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 899d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 900d5b485abSSatish Balay } 901d5b485abSSatish Balay 902d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 903d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 904d5b485abSSatish Balay j = pa[i]; 905d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 906d5b485abSSatish Balay msz += w1[j]; 907d5b485abSSatish Balay } 908d5b485abSSatish Balay 909c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 910a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 911c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 912d5b485abSSatish Balay 913c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 914c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 915c7dd2924SSatish Balay 916c7dd2924SSatish Balay ierr = PetscFree(onodes1);CHKERRQ(ierr); 917c7dd2924SSatish Balay ierr = PetscFree(olengths1);CHKERRQ(ierr); 918d5b485abSSatish Balay 919d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 920d5b485abSSatish Balay len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 92182502324SSatish Balay ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr); 922d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 923549d3d68SSatish Balay ierr = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr); 924d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 925d5b485abSSatish Balay tmp = (int*)(ptr + size); 926d5b485abSSatish Balay ctr = tmp + 2*msz; 927d5b485abSSatish Balay 928d5b485abSSatish Balay { 929d5b485abSSatish Balay int *iptr = tmp,ict = 0; 930d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 931d5b485abSSatish Balay j = pa[i]; 932d5b485abSSatish Balay iptr += ict; 933d5b485abSSatish Balay sbuf1[j] = iptr; 934d5b485abSSatish Balay ict = w1[j]; 935d5b485abSSatish Balay } 936d5b485abSSatish Balay } 937d5b485abSSatish Balay 938d5b485abSSatish Balay /* Form the outgoing messages */ 939d5b485abSSatish Balay /* Initialise the header space */ 940d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 941d5b485abSSatish Balay j = pa[i]; 942d5b485abSSatish Balay sbuf1[j][0] = 0; 943549d3d68SSatish Balay ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 944d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 945d5b485abSSatish Balay } 946d5b485abSSatish Balay 947d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 948d5b485abSSatish Balay for (i=0; i<ismax; i++) { 949549d3d68SSatish Balay ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 950d5b485abSSatish Balay irow_i = irow[i]; 951d5b485abSSatish Balay jmax = nrow[i]; 952d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 953d5b485abSSatish Balay row = irow_i[j]; 954233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 955233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 956233c2ff6SSatish Balay #else 957d5b485abSSatish Balay proc = rtable[row]; 958233c2ff6SSatish Balay #endif 959d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 960d5b485abSSatish Balay ctr[proc]++; 961d5b485abSSatish Balay *ptr[proc] = row; 962d5b485abSSatish Balay ptr[proc]++; 963d5b485abSSatish Balay } 964d5b485abSSatish Balay } 965d5b485abSSatish Balay /* Update the headers for the current IS */ 966d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 96706ef35abSSatish Balay if ((ctr_j = ctr[j])) { 968d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 969d5b485abSSatish Balay k = ++sbuf1_j[0]; 970d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 971d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 972d5b485abSSatish Balay } 973d5b485abSSatish Balay } 974d5b485abSSatish Balay } 975d5b485abSSatish Balay 976d5b485abSSatish Balay /* Now post the sends */ 977b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 978d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 979d5b485abSSatish Balay j = pa[i]; 980ca161407SBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 981d5b485abSSatish Balay } 982d5b485abSSatish Balay 983d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 984b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 98582502324SSatish Balay ierr = PetscMalloc((nrqs+1)*sizeof(int *),&rbuf2);CHKERRQ(ierr); 986d5b485abSSatish Balay rbuf2[0] = tmp + msz; 987d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 988d5b485abSSatish Balay j = pa[i]; 989d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 990d5b485abSSatish Balay } 991d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 992d5b485abSSatish Balay j = pa[i]; 993ca161407SBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 994d5b485abSSatish Balay } 995d5b485abSSatish Balay 996d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 997d5b485abSSatish Balay 998d5b485abSSatish Balay /* Receive messages*/ 999b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 1000b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 1001d5b485abSSatish Balay len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*); 100282502324SSatish Balay ierr = PetscMalloc(len,&sbuf2);CHKERRQ(ierr); 1003d5b485abSSatish Balay req_size = (int*)(sbuf2 + nrqr); 1004d5b485abSSatish Balay req_source = req_size + nrqr; 1005d5b485abSSatish Balay 1006d5b485abSSatish Balay { 1007df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 1008ca161407SBarry Smith int *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 1009d5b485abSSatish Balay 1010d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1011999d9058SBarry Smith ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 1012999d9058SBarry Smith req_size[idex] = 0; 1013999d9058SBarry Smith rbuf1_i = rbuf1[idex]; 1014d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 1015ca161407SBarry Smith ierr = MPI_Get_count(r_status1+i,MPI_INT,&end);CHKERRQ(ierr); 1016999d9058SBarry Smith ierr = PetscMalloc(end*sizeof(int),&sbuf2[idex]);CHKERRQ(ierr); 1017999d9058SBarry Smith sbuf2_i = sbuf2[idex]; 1018d5b485abSSatish Balay for (j=start; j<end; j++) { 1019d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 1020d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 1021d5b485abSSatish Balay sbuf2_i[j] = ncols; 1022999d9058SBarry Smith req_size[idex] += ncols; 1023d5b485abSSatish Balay } 1024999d9058SBarry Smith req_source[idex] = r_status1[i].MPI_SOURCE; 1025d5b485abSSatish Balay /* form the header */ 1026999d9058SBarry Smith sbuf2_i[0] = req_size[idex]; 1027d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 1028999d9058SBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 1029d5b485abSSatish Balay } 1030d5b485abSSatish Balay } 1031606d414cSSatish Balay ierr = PetscFree(r_status1);CHKERRQ(ierr); 1032606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 1033d5b485abSSatish Balay 1034d5b485abSSatish Balay /* recv buffer sizes */ 1035d5b485abSSatish Balay /* Receive messages*/ 1036d5b485abSSatish Balay 1037b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);CHKERRQ(ierr); 1038b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MatScalar*),&rbuf4);CHKERRQ(ierr); 1039b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 1040b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 1041b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 1042d5b485abSSatish Balay 1043d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1044999d9058SBarry Smith ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 1045999d9058SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*sizeof(int),&rbuf3[idex]);CHKERRQ(ierr); 1046999d9058SBarry Smith ierr = PetscMalloc(rbuf2[idex][0]*bs2*sizeof(MatScalar),&rbuf4[idex]);CHKERRQ(ierr); 1047999d9058SBarry Smith ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPI_INT, 1048999d9058SBarry Smith r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 1049999d9058SBarry Smith ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR, 1050999d9058SBarry Smith r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 1051d5b485abSSatish Balay } 1052606d414cSSatish Balay ierr = PetscFree(r_status2);CHKERRQ(ierr); 1053606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 1054d5b485abSSatish Balay 1055d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 1056b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 1057b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 1058d5b485abSSatish Balay 1059ca161407SBarry Smith ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 1060ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 1061606d414cSSatish Balay ierr = PetscFree(s_status1);CHKERRQ(ierr); 1062606d414cSSatish Balay ierr = PetscFree(s_status2);CHKERRQ(ierr); 1063606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 1064606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 1065d5b485abSSatish Balay 1066d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 106782502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(int *),&sbuf_aj);CHKERRQ(ierr); 1068d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 106982502324SSatish Balay ierr = PetscMalloc((j+1)*sizeof(int),&sbuf_aj[0]);CHKERRQ(ierr); 1070d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1071d5b485abSSatish Balay 1072b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 1073d5b485abSSatish Balay { 1074d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1075d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1076d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 1077d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 1078d5b485abSSatish Balay ct2 = 0; 1079d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1080d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 1081d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1082e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1083d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1084d5b485abSSatish Balay ncols = nzA + nzB; 1085d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1086d5b485abSSatish Balay 1087d5b485abSSatish Balay /* load the column indices for this row into cols*/ 1088d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 1089d5b485abSSatish Balay for (l=0; l<nzB; l++) { 1090d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1091d5b485abSSatish Balay else break; 1092d5b485abSSatish Balay } 1093d5b485abSSatish Balay imark = l; 1094d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1095d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1096d5b485abSSatish Balay ct2 += ncols; 1097d5b485abSSatish Balay } 1098d5b485abSSatish Balay } 1099ca161407SBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 1100d5b485abSSatish Balay } 1101d5b485abSSatish Balay } 1102b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 1103b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 1104d5b485abSSatish Balay 1105d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 110682502324SSatish Balay ierr = PetscMalloc((nrqr+1)*sizeof(MatScalar *),&sbuf_aa);CHKERRQ(ierr); 1107d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 110882502324SSatish Balay ierr = PetscMalloc((j+1)*bs2*sizeof(MatScalar),&sbuf_aa[0]);CHKERRQ(ierr); 1109a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1110d5b485abSSatish Balay 1111b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 1112d5b485abSSatish Balay { 1113d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1114d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1115d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 1116d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 1117d5b485abSSatish Balay ct2 = 0; 1118d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1119d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 1120d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1121e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1122d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1123d5b485abSSatish Balay ncols = nzA + nzB; 1124d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1125a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 1126d5b485abSSatish Balay 1127d5b485abSSatish Balay /* load the column values for this row into vals*/ 11285b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 1129d5b485abSSatish Balay for (l=0; l<nzB; l++) { 11303a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 11313eda8832SBarry Smith ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 11323a40ed3dSBarry Smith } 1133d5b485abSSatish Balay else break; 1134d5b485abSSatish Balay } 1135d5b485abSSatish Balay imark = l; 11363a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 11373eda8832SBarry Smith ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 11383a40ed3dSBarry Smith } 11393a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 11403eda8832SBarry Smith ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 11413a40ed3dSBarry Smith } 1142d5b485abSSatish Balay ct2 += ncols; 1143d5b485abSSatish Balay } 1144d5b485abSSatish Balay } 11453eda8832SBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 1146d5b485abSSatish Balay } 1147d5b485abSSatish Balay } 1148b0a32e0cSBarry Smith ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 1149b0a32e0cSBarry Smith ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 1150606d414cSSatish Balay ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1151d5b485abSSatish Balay 1152d5b485abSSatish Balay /* Form the matrix */ 1153d5b485abSSatish Balay /* create col map */ 1154d5b485abSSatish Balay { 1155d5b485abSSatish Balay int *icol_i; 1156233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1157233c2ff6SSatish Balay /* Create row map*/ 1158b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&colmaps);CHKERRQ(ierr); 1159ff0794f7SSatish Balay for (i=0; i<ismax; i++) { 1160ff0794f7SSatish Balay ierr = PetscTableCreate(ncol[i]+1,&colmaps[i]);CHKERRQ(ierr); 1161233c2ff6SSatish Balay } 1162233c2ff6SSatish Balay #else 1163a2feddc5SSatish Balay len = (1+ismax)*sizeof(int*)+ ismax*c->Nbs*sizeof(int); 116482502324SSatish Balay ierr = PetscMalloc(len,&cmap);CHKERRQ(ierr); 1165d5b485abSSatish Balay cmap[0] = (int *)(cmap + ismax); 1166549d3d68SSatish Balay ierr = PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int));CHKERRQ(ierr); 1167a2feddc5SSatish Balay for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; } 1168233c2ff6SSatish Balay #endif 1169d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1170d5b485abSSatish Balay jmax = ncol[i]; 1171d5b485abSSatish Balay icol_i = icol[i]; 1172233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1173233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1174233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1175001aedefSBarry Smith ierr = PetscTableAdd(lcol1_gcol1,icol_i[j]+1,j+1);CHKERRQ(ierr); 1176233c2ff6SSatish Balay } 1177233c2ff6SSatish Balay #else 1178d5b485abSSatish Balay cmap_i = cmap[i]; 1179d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1180d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1181d5b485abSSatish Balay } 1182233c2ff6SSatish Balay #endif 1183d5b485abSSatish Balay } 1184d5b485abSSatish Balay } 1185d5b485abSSatish Balay 1186d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1187d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1188d5b485abSSatish Balay len = (1+ismax)*sizeof(int*)+ j*sizeof(int); 118982502324SSatish Balay ierr = PetscMalloc(len,&lens);CHKERRQ(ierr); 1190d5b485abSSatish Balay lens[0] = (int *)(lens + ismax); 1191549d3d68SSatish Balay ierr = PetscMemzero(lens[0],j*sizeof(int));CHKERRQ(ierr); 1192d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1193d5b485abSSatish Balay 1194d5b485abSSatish Balay /* Update lens from local data */ 1195d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1196d5b485abSSatish Balay jmax = nrow[i]; 1197233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1198233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1199233c2ff6SSatish Balay #else 1200d5b485abSSatish Balay cmap_i = cmap[i]; 1201233c2ff6SSatish Balay #endif 1202d5b485abSSatish Balay irow_i = irow[i]; 1203d5b485abSSatish Balay lens_i = lens[i]; 1204d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1205d5b485abSSatish Balay row = irow_i[j]; 1206233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1207233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1208233c2ff6SSatish Balay #else 1209d5b485abSSatish Balay proc = rtable[row]; 1210233c2ff6SSatish Balay #endif 1211d5b485abSSatish Balay if (proc == rank) { 121236f4e84dSSatish Balay /* Get indices from matA and then from matB */ 121336f4e84dSSatish Balay row = row - rstart; 121436f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 121536f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1216233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1217233c2ff6SSatish Balay for (k=0; k<nzA; k++) { 1218233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1219233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1220233c2ff6SSatish Balay } 1221233c2ff6SSatish Balay for (k=0; k<nzB; k++) { 1222233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1223233c2ff6SSatish Balay if (tt) { lens_i[j]++; } 1224233c2ff6SSatish Balay } 1225233c2ff6SSatish Balay #else 1226ca161407SBarry Smith for (k=0; k<nzA; k++) { 122736f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++; } 1228ca161407SBarry Smith } 1229ca161407SBarry Smith for (k=0; k<nzB; k++) { 123036f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++; } 1231d5b485abSSatish Balay } 1232233c2ff6SSatish Balay #endif 1233d5b485abSSatish Balay } 1234a2feddc5SSatish Balay } 1235ca161407SBarry Smith } 1236233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1237d5b485abSSatish Balay /* Create row map*/ 1238b0a32e0cSBarry Smith ierr = PetscMalloc((1+ismax)*sizeof(PetscTable),&rowmaps);CHKERRQ(ierr); 1239ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 1240ff0794f7SSatish Balay ierr = PetscTableCreate(nrow[i]+1,&rowmaps[i]);CHKERRQ(ierr); 1241233c2ff6SSatish Balay } 1242233c2ff6SSatish Balay #else 1243233c2ff6SSatish Balay /* Create row map*/ 1244233c2ff6SSatish Balay len = (1+ismax)*sizeof(int*)+ ismax*Mbs*sizeof(int); 124582502324SSatish Balay ierr = PetscMalloc(len,&rmap);CHKERRQ(ierr); 1246d5b485abSSatish Balay rmap[0] = (int *)(rmap + ismax); 1247233c2ff6SSatish Balay ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(int));CHKERRQ(ierr); 1248233c2ff6SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + Mbs;} 1249233c2ff6SSatish Balay #endif 1250d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1251d5b485abSSatish Balay irow_i = irow[i]; 1252d5b485abSSatish Balay jmax = nrow[i]; 1253233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1254233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1255233c2ff6SSatish Balay for (j=0; j<jmax; j++) { 1256233c2ff6SSatish Balay ierr = PetscTableAdd(lrow1_grow1,irow_i[j]+1,j+1);CHKERRQ(ierr); 1257233c2ff6SSatish Balay } 1258233c2ff6SSatish Balay #else 1259233c2ff6SSatish Balay rmap_i = rmap[i]; 1260d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1261d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1262d5b485abSSatish Balay } 1263233c2ff6SSatish Balay #endif 1264d5b485abSSatish Balay } 1265d5b485abSSatish Balay 1266d5b485abSSatish Balay /* Update lens from offproc data */ 1267d5b485abSSatish Balay { 1268d5b485abSSatish Balay int *rbuf2_i,*rbuf3_i,*sbuf1_i; 1269d5b485abSSatish Balay 1270d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1271ca161407SBarry Smith ierr = MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);CHKERRQ(ierr); 1272999d9058SBarry Smith idex = pa[i]; 1273999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1274d5b485abSSatish Balay jmax = sbuf1_i[0]; 1275d5b485abSSatish Balay ct1 = 2*jmax+1; 1276d5b485abSSatish Balay ct2 = 0; 1277d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1278d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1279d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1280d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1281d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1282d5b485abSSatish Balay lens_i = lens[is_no]; 1283233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1284233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1285233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1286233c2ff6SSatish Balay #else 1287d5b485abSSatish Balay cmap_i = cmap[is_no]; 1288d5b485abSSatish Balay rmap_i = rmap[is_no]; 1289233c2ff6SSatish Balay #endif 1290d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1291233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1292233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1293233c2ff6SSatish Balay row--; 129429bbc08cSBarry Smith if(row < 0) { SETERRQ(1,"row not found in table"); } 1295233c2ff6SSatish Balay #else 1296d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1297233c2ff6SSatish Balay #endif 1298d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1299d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1300233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1301233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 1302233c2ff6SSatish Balay if (tt) { 1303233c2ff6SSatish Balay lens_i[row]++; 1304233c2ff6SSatish Balay } 1305233c2ff6SSatish Balay #else 1306d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1307d5b485abSSatish Balay lens_i[row]++; 1308d5b485abSSatish Balay } 1309233c2ff6SSatish Balay #endif 1310d5b485abSSatish Balay } 1311d5b485abSSatish Balay } 1312d5b485abSSatish Balay } 1313d5b485abSSatish Balay } 1314d5b485abSSatish Balay } 1315606d414cSSatish Balay ierr = PetscFree(r_status3);CHKERRQ(ierr); 1316606d414cSSatish Balay ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1317ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr); 1318606d414cSSatish Balay ierr = PetscFree(s_status3);CHKERRQ(ierr); 1319606d414cSSatish Balay ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1320d5b485abSSatish Balay 1321d5b485abSSatish Balay /* Create the submatrices */ 1322d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1323d5b485abSSatish Balay /* 1324d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1325d5b485abSSatish Balay */ 1326d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1327df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 132836f4e84dSSatish Balay if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) { 132929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1330d5b485abSSatish Balay } 13310f5bd95cSBarry Smith ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(int),&flag);CHKERRQ(ierr); 13320f5bd95cSBarry Smith if (flag == PETSC_FALSE) { 133329bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 1334d5b485abSSatish Balay } 1335d5b485abSSatish Balay /* Initial matrix as if empty */ 1336549d3d68SSatish Balay ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(int));CHKERRQ(ierr); 1337ce60de66SLois Curfman McInnes submats[i]->factor = C->factor; 1338d5b485abSSatish Balay } 1339ca161407SBarry Smith } else { 1340d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1341e7ef3d9dSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,nrow[i]*bs,ncol[i]*bs,nrow[i]*bs,ncol[i]*bs,submats+i);CHKERRQ(ierr); 1342e7ef3d9dSBarry Smith ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr); 1343e7ef3d9dSBarry Smith ierr = MatSeqBAIJSetPreallocation(submats[i],a->bs,0,lens[i]);CHKERRQ(ierr); 1344e7ef3d9dSBarry Smith ierr = MatSeqSBAIJSetPreallocation(submats[i],a->bs,0,lens[i]);CHKERRQ(ierr); 1345d5b485abSSatish Balay } 1346d5b485abSSatish Balay } 1347d5b485abSSatish Balay 1348d5b485abSSatish Balay /* Assemble the matrices */ 1349d5b485abSSatish Balay /* First assemble the local rows */ 1350d5b485abSSatish Balay { 135136f4e84dSSatish Balay int ilen_row,*imat_ilen,*imat_j,*imat_i; 13523eda8832SBarry Smith MatScalar *imat_a; 1353d5b485abSSatish Balay 1354d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1355df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[i]->data; 1356d5b485abSSatish Balay imat_ilen = mat->ilen; 1357d5b485abSSatish Balay imat_j = mat->j; 1358d5b485abSSatish Balay imat_i = mat->i; 1359d5b485abSSatish Balay imat_a = mat->a; 1360233c2ff6SSatish Balay 1361233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1362233c2ff6SSatish Balay lcol1_gcol1 = colmaps[i]; 1363233c2ff6SSatish Balay lrow1_grow1 = rowmaps[i]; 1364233c2ff6SSatish Balay #else 1365d5b485abSSatish Balay cmap_i = cmap[i]; 1366d5b485abSSatish Balay rmap_i = rmap[i]; 1367233c2ff6SSatish Balay #endif 1368d5b485abSSatish Balay irow_i = irow[i]; 1369d5b485abSSatish Balay jmax = nrow[i]; 1370d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1371d5b485abSSatish Balay row = irow_i[j]; 1372233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1373233c2ff6SSatish Balay ierr = PetscGetProc(row,size,c->rowners,&proc);CHKERRQ(ierr); 1374233c2ff6SSatish Balay #else 1375d5b485abSSatish Balay proc = rtable[row]; 1376233c2ff6SSatish Balay #endif 1377d5b485abSSatish Balay if (proc == rank) { 137836f4e84dSSatish Balay row = row - rstart; 137936f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 138036f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 138136f4e84dSSatish Balay cworkA = a_j + a_i[row]; 138236f4e84dSSatish Balay cworkB = b_j + b_i[row]; 138336f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 138436f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 1385233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1386233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+rstart+1,&row);CHKERRQ(ierr); 1387233c2ff6SSatish Balay row--; 138829bbc08cSBarry Smith if (row < 0) { SETERRQ(1,"row not found in table"); } 1389233c2ff6SSatish Balay #else 139036f4e84dSSatish Balay row = rmap_i[row + rstart]; 1391233c2ff6SSatish Balay #endif 1392df36731bSBarry Smith mat_i = imat_i[row]; 139336f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1394d5b485abSSatish Balay mat_j = imat_j + mat_i; 139536f4e84dSSatish Balay ilen_row = imat_ilen[row]; 139636f4e84dSSatish Balay 139736f4e84dSSatish Balay /* load the column indices for this row into cols*/ 139836f4e84dSSatish Balay for (l=0; l<nzB; l++) { 139936f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 1400233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1401233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,ctmp+1,&tcol);CHKERRQ(ierr); 1402233c2ff6SSatish Balay if (tcol) { 1403233c2ff6SSatish Balay #else 140436f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1405233c2ff6SSatish Balay #endif 1406df36731bSBarry Smith *mat_j++ = tcol - 1; 14073eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1408549d3d68SSatish Balay mat_a += bs2; 1409d5b485abSSatish Balay ilen_row++; 1410d5b485abSSatish Balay } 1411ca161407SBarry Smith } else break; 141236f4e84dSSatish Balay } 141336f4e84dSSatish Balay imark = l; 141436f4e84dSSatish Balay for (l=0; l<nzA; l++) { 1415233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1416233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1417233c2ff6SSatish Balay if (tcol) { 1418233c2ff6SSatish Balay #else 141936f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 1420233c2ff6SSatish Balay #endif 142136f4e84dSSatish Balay *mat_j++ = tcol - 1; 14223eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1423549d3d68SSatish Balay mat_a += bs2; 142436f4e84dSSatish Balay ilen_row++; 142536f4e84dSSatish Balay } 142636f4e84dSSatish Balay } 142736f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 1428233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1429233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1430233c2ff6SSatish Balay if (tcol) { 1431233c2ff6SSatish Balay #else 143236f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1433233c2ff6SSatish Balay #endif 143436f4e84dSSatish Balay *mat_j++ = tcol - 1; 14353eda8832SBarry Smith ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1436549d3d68SSatish Balay mat_a += bs2; 143736f4e84dSSatish Balay ilen_row++; 143836f4e84dSSatish Balay } 143936f4e84dSSatish Balay } 1440d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1441d5b485abSSatish Balay } 1442d5b485abSSatish Balay } 144336f4e84dSSatish Balay 1444d5b485abSSatish Balay } 1445d5b485abSSatish Balay } 1446d5b485abSSatish Balay 1447d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1448d5b485abSSatish Balay { 1449d5b485abSSatish Balay int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1450d5b485abSSatish Balay int *imat_j,*imat_i; 14513eda8832SBarry Smith MatScalar *imat_a,*rbuf4_i; 1452d5b485abSSatish Balay 1453d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1454ca161407SBarry Smith ierr = MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);CHKERRQ(ierr); 1455999d9058SBarry Smith idex = pa[i]; 1456999d9058SBarry Smith sbuf1_i = sbuf1[idex]; 1457d5b485abSSatish Balay jmax = sbuf1_i[0]; 1458d5b485abSSatish Balay ct1 = 2*jmax + 1; 1459d5b485abSSatish Balay ct2 = 0; 1460d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1461d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1462d5b485abSSatish Balay rbuf4_i = rbuf4[i]; 1463d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1464d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1465233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1466233c2ff6SSatish Balay lrow1_grow1 = rowmaps[is_no]; 1467233c2ff6SSatish Balay lcol1_gcol1 = colmaps[is_no]; 1468233c2ff6SSatish Balay #else 1469d5b485abSSatish Balay rmap_i = rmap[is_no]; 1470d5b485abSSatish Balay cmap_i = cmap[is_no]; 1471233c2ff6SSatish Balay #endif 1472df36731bSBarry Smith mat = (Mat_SeqBAIJ*)submats[is_no]->data; 1473d5b485abSSatish Balay imat_ilen = mat->ilen; 1474d5b485abSSatish Balay imat_j = mat->j; 1475d5b485abSSatish Balay imat_i = mat->i; 1476d5b485abSSatish Balay imat_a = mat->a; 1477d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1478d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1479d5b485abSSatish Balay row = sbuf1_i[ct1]; 1480233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1481233c2ff6SSatish Balay ierr = PetscTableFind(lrow1_grow1,row+1,&row);CHKERRQ(ierr); 1482233c2ff6SSatish Balay row--; 148329bbc08cSBarry Smith if(row < 0) { SETERRQ(1,"row not found in table"); } 1484233c2ff6SSatish Balay #else 1485d5b485abSSatish Balay row = rmap_i[row]; 1486233c2ff6SSatish Balay #endif 1487d5b485abSSatish Balay ilen = imat_ilen[row]; 1488df36731bSBarry Smith mat_i = imat_i[row]; 148936f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1490d5b485abSSatish Balay mat_j = imat_j + mat_i; 1491d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1492d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1493233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1494233c2ff6SSatish Balay ierr = PetscTableFind(lcol1_gcol1,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1495233c2ff6SSatish Balay if (tcol) { 1496233c2ff6SSatish Balay #else 1497d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1498233c2ff6SSatish Balay #endif 1499df36731bSBarry Smith *mat_j++ = tcol - 1; 150036f4e84dSSatish Balay /* *mat_a++= rbuf4_i[ct2]; */ 15013eda8832SBarry Smith ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1502549d3d68SSatish Balay mat_a += bs2; 1503d5b485abSSatish Balay ilen++; 1504d5b485abSSatish Balay } 1505d5b485abSSatish Balay } 1506d5b485abSSatish Balay imat_ilen[row] = ilen; 1507d5b485abSSatish Balay } 1508d5b485abSSatish Balay } 1509d5b485abSSatish Balay } 1510d5b485abSSatish Balay } 1511606d414cSSatish Balay ierr = PetscFree(r_status4);CHKERRQ(ierr); 1512606d414cSSatish Balay ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1513ca161407SBarry Smith ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr); 1514606d414cSSatish Balay ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1515606d414cSSatish Balay ierr = PetscFree(s_status4);CHKERRQ(ierr); 1516d5b485abSSatish Balay 1517d5b485abSSatish Balay /* Restore the indices */ 1518d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1519d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1520d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1521d5b485abSSatish Balay } 1522d5b485abSSatish Balay 1523d5b485abSSatish Balay /* Destroy allocated memory */ 1524606d414cSSatish Balay ierr = PetscFree(irow);CHKERRQ(ierr); 1525606d414cSSatish Balay ierr = PetscFree(w1);CHKERRQ(ierr); 1526606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 1527d5b485abSSatish Balay 1528606d414cSSatish Balay ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1529606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1530d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1531606d414cSSatish Balay ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1532d5b485abSSatish Balay } 1533d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1534606d414cSSatish Balay ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1535606d414cSSatish Balay ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1536d5b485abSSatish Balay } 1537d5b485abSSatish Balay 1538606d414cSSatish Balay ierr = PetscFree(sbuf2);CHKERRQ(ierr); 1539606d414cSSatish Balay ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1540606d414cSSatish Balay ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1541606d414cSSatish Balay ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1542606d414cSSatish Balay ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1543606d414cSSatish Balay ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1544606d414cSSatish Balay ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1545d5b485abSSatish Balay 1546233c2ff6SSatish Balay #if defined (PETSC_USE_CTABLE) 1547ff0794f7SSatish Balay for (i=0; i<ismax; i++){ 1548233c2ff6SSatish Balay ierr = PetscTableDelete(rowmaps[i]);CHKERRQ(ierr); 1549233c2ff6SSatish Balay ierr = PetscTableDelete(colmaps[i]);CHKERRQ(ierr); 1550233c2ff6SSatish Balay } 1551233c2ff6SSatish Balay ierr = PetscFree(colmaps);CHKERRQ(ierr); 1552233c2ff6SSatish Balay ierr = PetscFree(rowmaps);CHKERRQ(ierr); 1553233c2ff6SSatish Balay #else 1554606d414cSSatish Balay ierr = PetscFree(rmap);CHKERRQ(ierr); 1555233c2ff6SSatish Balay ierr = PetscFree(cmap);CHKERRQ(ierr); 1556233c2ff6SSatish Balay #endif 1557606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1558d5b485abSSatish Balay 1559d5b485abSSatish Balay for (i=0; i<ismax; i++) { 156036f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 156136f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1562d5b485abSSatish Balay } 156334e6ae18SSatish Balay 15643a40ed3dSBarry Smith PetscFunctionReturn(0); 1565d5b485abSSatish Balay } 1566ca161407SBarry Smith 1567