1d5b485abSSatish Balay #ifndef lint 2*70f55243SBarry Smith static char vcid[] = "$Id: baijov.c,v 1.9 1996/08/06 14:36:31 balay Exp bsmith $"; 3d5b485abSSatish Balay #endif 4d5b485abSSatish Balay /* 5d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 6d5b485abSSatish Balay and to find submatrices that were shared across processors. 7d5b485abSSatish Balay */ 8*70f55243SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h" 9d5b485abSSatish Balay #include "src/inline/bitarray.h" 10d5b485abSSatish Balay 11df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat, int, IS *); 12df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat , int , char **,int*, int**); 13df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat , int, int **, int**, int* ); 14df36731bSBarry Smith extern int MatGetRow_MPIBAIJ(Mat,int,int*,int**,Scalar**); 15df36731bSBarry Smith extern int MatRestoreRow_MPIBAIJ(Mat,int,int*,int**,Scalar**); 16d5b485abSSatish Balay 17e757655aSSatish Balay 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; 21df36731bSBarry Smith int ierr,isz,bs = baij->bs,Nbs,n,i,j,*idx,*nidx,ival; 22df36731bSBarry Smith char *table; 23df36731bSBarry Smith 24df36731bSBarry Smith Nbs = baij->Nbs; 25df36731bSBarry Smith table = (char *) PetscMalloc((Nbs/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 26df36731bSBarry Smith nidx = (int *) PetscMalloc((Nbs+1)*sizeof(int)); CHKPTRQ(nidx); 27df36731bSBarry Smith 28df36731bSBarry Smith for (i=0; i<imax; i++) { 29df36731bSBarry Smith isz = 0; 30df36731bSBarry Smith PetscMemzero(table,(Nbs/BITSPERBYTE +1)*sizeof(char)); 3136f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 3236f4e84dSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 33e757655aSSatish Balay for (j=0; j<n ; j++) { 34df36731bSBarry Smith ival = idx[j]/bs; /* convert the indices into block indices */ 35e757655aSSatish Balay if (ival>Nbs) SETERRQ(1,"MatCompressIndicesGeneral_MPIBAIJ: index greater than mat-dim"); 36df36731bSBarry Smith if(!BT_LOOKUP(table, ival)) { nidx[isz++] = ival;} 37df36731bSBarry Smith } 3836f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 3936f4e84dSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is_out+i)); CHKERRQ(ierr); 40df36731bSBarry Smith } 41df36731bSBarry Smith PetscFree(table); 42df36731bSBarry Smith PetscFree(nidx); 43df36731bSBarry Smith return 0; 44df36731bSBarry Smith } 45df36731bSBarry Smith 46e757655aSSatish Balay static int MatCompressIndicesSorted_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 47e757655aSSatish Balay { 48e757655aSSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 49e757655aSSatish Balay int ierr,bs=baij->bs,i,j,k,val,n,*idx,*nidx,Nbs=baij->Nbs,*idx_local; 50e757655aSSatish Balay PetscTruth flg; 51e757655aSSatish Balay 52e757655aSSatish Balay for (i=0; i<imax; i++) { 53e757655aSSatish Balay ierr = ISSorted(is_in[i],&flg); CHKERRQ(ierr); 54e757655aSSatish Balay if (!flg) SETERRQ(1,"MatCompressIndicesSorted_MPIBAIJ: Indices are not sorted"); 55e757655aSSatish Balay } 56e757655aSSatish Balay nidx = (int *) PetscMalloc((Nbs+1)*sizeof(int)); CHKPTRQ(nidx); 57e757655aSSatish Balay /* Now chech if the indices are in block order */ 58e757655aSSatish Balay for (i=0; i<imax; i++) { 59e757655aSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 60e757655aSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 61e757655aSSatish Balay if ( n%bs !=0 ) SETERRA(1,"MatCompressIndicesSorted_MPIBAIJ: Indices are not block ordered"); 62e757655aSSatish Balay 63e757655aSSatish Balay n = n/bs; /* The reduced index size */ 64e757655aSSatish Balay idx_local = idx; 65e757655aSSatish Balay for (j=0; j<n ; j++) { 66e757655aSSatish Balay val = idx_local[0]; 67e757655aSSatish Balay if(val%bs != 0) SETERRA(1,"MatCompressIncicesSorted_MPIBAIJ: Indices are not block ordered"); 68e757655aSSatish Balay for (k=0; k<bs; k++) { 69e757655aSSatish Balay if ( val+k != idx_local[k]) SETERRA(1,"MatCompressIncicesSorted_MPIBAIJ: Indices are not block ordered"); 70e757655aSSatish Balay } 71e757655aSSatish Balay nidx[j] = val/bs; 72e757655aSSatish Balay idx_local +=bs; 73e757655aSSatish Balay } 74e757655aSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 75e757655aSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF,n,nidx,(is_out+i)); CHKERRQ(ierr); 76e757655aSSatish Balay } 77e757655aSSatish Balay PetscFree(nidx); 78e757655aSSatish Balay return 0; 79e757655aSSatish Balay } 80e757655aSSatish Balay 8136f4e84dSSatish Balay static int MatExpandIndices_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 82df36731bSBarry Smith { 83df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 84df36731bSBarry Smith int ierr,bs = baij->bs,Nbs,n,i,j,k,*idx,*nidx; 85df36731bSBarry Smith 86df36731bSBarry Smith Nbs = baij->Nbs; 87df36731bSBarry Smith 88df36731bSBarry Smith nidx = (int *) PetscMalloc((Nbs*bs+1)*sizeof(int)); CHKPTRQ(nidx); 89df36731bSBarry Smith 90df36731bSBarry Smith for ( i=0; i<imax; i++ ) { 9136f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 9236f4e84dSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 93df36731bSBarry Smith for (j=0; j<n ; ++j){ 94df36731bSBarry Smith for (k=0; k<bs; k++) 95df36731bSBarry Smith nidx[j*bs+k] = idx[j]*bs+k; 96df36731bSBarry Smith } 9736f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 9836f4e84dSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, n*bs, nidx, (is_out+i)); CHKERRQ(ierr); 99df36731bSBarry Smith } 100df36731bSBarry Smith PetscFree(nidx); 101df36731bSBarry Smith return 0; 102df36731bSBarry Smith } 103df36731bSBarry Smith 104df36731bSBarry Smith 105df36731bSBarry Smith int MatIncreaseOverlap_MPIBAIJ(Mat C, int imax, IS *is, int ov) 106d5b485abSSatish Balay { 107d5b485abSSatish Balay int i, ierr; 10836f4e84dSSatish Balay IS *is_new; 10936f4e84dSSatish Balay 11036f4e84dSSatish Balay is_new = (IS *)PetscMalloc(imax*sizeof(IS)); CHKPTRQ(is_new); 111df36731bSBarry Smith /* Convert the indices into block format */ 112e757655aSSatish Balay ierr = MatCompressIndicesGeneral_MPIBAIJ(C, imax, is,is_new); CHKERRQ(ierr); 113df36731bSBarry Smith if (ov < 0){ SETERRQ(1," MatIncreaseOverlap_MPIBAIJ:Negative overlap specified\n");} 114d5b485abSSatish Balay for (i=0; i<ov; ++i) { 11536f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new); CHKERRQ(ierr); 116d5b485abSSatish Balay } 11736f4e84dSSatish Balay for (i=0; i<imax; i++) ISDestroy(is[i]); 11836f4e84dSSatish Balay ierr = MatExpandIndices_MPIBAIJ(C, imax, is_new,is); CHKERRQ(ierr); 11936f4e84dSSatish Balay for (i=0; i<imax; i++) ISDestroy(is_new[i]); 12036f4e84dSSatish Balay PetscFree(is_new); 121d5b485abSSatish Balay return 0; 122d5b485abSSatish Balay } 123d5b485abSSatish Balay 124d5b485abSSatish Balay /* 125d5b485abSSatish Balay Sample message format: 126d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 127d5b485abSSatish Balay to index sets 1s[1], is[5] 128d5b485abSSatish Balay mesg [0] = 2 ( no of index sets in the mesg) 129d5b485abSSatish Balay ----------- 130d5b485abSSatish Balay mesg [1] = 1 => is[1] 131d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 132d5b485abSSatish Balay ----------- 133d5b485abSSatish Balay mesg [5] = 5 => is[5] 134d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 135d5b485abSSatish Balay ----------- 136d5b485abSSatish Balay mesg [7] 137d5b485abSSatish Balay mesg [n] datas[1] 138d5b485abSSatish Balay ----------- 139d5b485abSSatish Balay mesg[n+1] 140d5b485abSSatish Balay mesg[m] data(is[5]) 141d5b485abSSatish Balay ----------- 142d5b485abSSatish Balay 143d5b485abSSatish Balay Notes: 144d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 145d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 146d5b485abSSatish Balay */ 147df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat C, int imax, IS *is) 148d5b485abSSatish Balay { 149df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 150d5b485abSSatish Balay int **idx, *n, *w1, *w2, *w3, *w4, *rtable,**data,len,*idx_i; 151df36731bSBarry Smith int size,rank,Mbs,i,j,k,ierr,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 152d5b485abSSatish Balay int *ctr,*pa,tag,*tmp,bsz,nrqr,*isz,*isz1,**xdata,bsz1,**rbuf2; 153d5b485abSSatish Balay char **table; 154d5b485abSSatish Balay MPI_Comm comm; 155d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 156d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 157d5b485abSSatish Balay 158d5b485abSSatish Balay comm = C->comm; 159d5b485abSSatish Balay tag = C->tag; 160d5b485abSSatish Balay size = c->size; 161d5b485abSSatish Balay rank = c->rank; 162df36731bSBarry Smith Mbs = c->Mbs; 163d5b485abSSatish Balay 164df36731bSBarry Smith len = (imax+1)*sizeof(int *) + (imax + Mbs)*sizeof(int); 165d5b485abSSatish Balay idx = (int **) PetscMalloc(len); CHKPTRQ(idx); 166d5b485abSSatish Balay n = (int *) (idx + imax); 167d5b485abSSatish Balay rtable = n + imax; 168d5b485abSSatish Balay 169d5b485abSSatish Balay for (i=0; i<imax; i++) { 170d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]); CHKERRQ(ierr); 171d5b485abSSatish Balay ierr = ISGetSize(is[i],&n[i]); CHKERRQ(ierr); 172d5b485abSSatish Balay } 173d5b485abSSatish Balay 174d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 175d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 176d5b485abSSatish Balay len = c->rowners[i+1]; 177d5b485abSSatish Balay for (; j<len; j++) { 178d5b485abSSatish Balay rtable[j] = i; 179d5b485abSSatish Balay } 180d5b485abSSatish Balay } 181d5b485abSSatish Balay 182d5b485abSSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 183d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 184d5b485abSSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1);/* mesg size */ 185d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 186d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 187d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 188d5b485abSSatish Balay PetscMemzero(w1,size*3*sizeof(int)); /* initialise work vector*/ 189d5b485abSSatish Balay for (i=0; i<imax; i++) { 190d5b485abSSatish Balay PetscMemzero(w4,size*sizeof(int)); /* initialise work vector*/ 191d5b485abSSatish Balay idx_i = idx[i]; 192d5b485abSSatish Balay len = n[i]; 193d5b485abSSatish Balay for (j=0; j<len; j++) { 194d5b485abSSatish Balay row = idx_i[j]; 195d5b485abSSatish Balay proc = rtable[row]; 196d5b485abSSatish Balay w4[proc]++; 197d5b485abSSatish Balay } 198d5b485abSSatish Balay for (j=0; j<size; j++){ 199d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 200d5b485abSSatish Balay } 201d5b485abSSatish Balay } 202d5b485abSSatish Balay 203d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 204d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 205d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 206d5b485abSSatish Balay w3[rank] = 0; 207d5b485abSSatish Balay for ( i=0; i<size; i++) { 208d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 209d5b485abSSatish Balay } 210d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 211d5b485abSSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); 212d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 213d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 214d5b485abSSatish Balay } 215d5b485abSSatish Balay 216d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 217d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 218d5b485abSSatish Balay j = pa[i]; 219d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 220d5b485abSSatish Balay msz += w1[j]; 221d5b485abSSatish Balay } 222d5b485abSSatish Balay 223d5b485abSSatish Balay 224d5b485abSSatish Balay /* Do a global reduction to determine how many messages to expect*/ 225d5b485abSSatish Balay { 226d5b485abSSatish Balay int *rw1, *rw2; 227d5b485abSSatish Balay rw1 = (int *) PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 228d5b485abSSatish Balay rw2 = rw1+size; 229d5b485abSSatish Balay MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm); 230d5b485abSSatish Balay bsz = rw1[rank]; 231d5b485abSSatish Balay MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm); 232d5b485abSSatish Balay nrqr = rw2[rank]; 233d5b485abSSatish Balay PetscFree(rw1); 234d5b485abSSatish Balay } 235d5b485abSSatish Balay 236d5b485abSSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 237d5b485abSSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 238d5b485abSSatish Balay rbuf = (int**) PetscMalloc(len); CHKPTRQ(rbuf); 239d5b485abSSatish Balay rbuf[0] = (int *) (rbuf + nrqr); 240d5b485abSSatish Balay for (i=1; i<nrqr; ++i) rbuf[i] = rbuf[i-1] + bsz; 241d5b485abSSatish Balay 242d5b485abSSatish Balay /* Post the receives */ 243d5b485abSSatish Balay r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 244d5b485abSSatish Balay CHKPTRQ(r_waits1); 245d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 246d5b485abSSatish Balay MPI_Irecv(rbuf[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i); 247d5b485abSSatish Balay } 248d5b485abSSatish Balay 249d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 250d5b485abSSatish Balay len = 2*size*sizeof(int*) + (size+msz)*sizeof(int); 251d5b485abSSatish Balay outdat = (int **)PetscMalloc(len); CHKPTRQ(outdat); 252d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 253d5b485abSSatish Balay PetscMemzero(outdat,2*size*sizeof(int*)); 254d5b485abSSatish Balay tmp = (int *) (outdat + 2*size); 255d5b485abSSatish Balay ctr = tmp + msz; 256d5b485abSSatish Balay 257d5b485abSSatish Balay { 258d5b485abSSatish Balay int *iptr = tmp,ict = 0; 259d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 260d5b485abSSatish Balay j = pa[i]; 261d5b485abSSatish Balay iptr += ict; 262d5b485abSSatish Balay outdat[j] = iptr; 263d5b485abSSatish Balay ict = w1[j]; 264d5b485abSSatish Balay } 265d5b485abSSatish Balay } 266d5b485abSSatish Balay 267d5b485abSSatish Balay /* Form the outgoing messages */ 268d5b485abSSatish Balay /*plug in the headers*/ 269d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 270d5b485abSSatish Balay j = pa[i]; 271d5b485abSSatish Balay outdat[j][0] = 0; 272d5b485abSSatish Balay PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int)); 273d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 274d5b485abSSatish Balay } 275d5b485abSSatish Balay 276d5b485abSSatish Balay /* Memory for doing local proc's work*/ 277d5b485abSSatish Balay { 278d5b485abSSatish Balay int *d_p; 279d5b485abSSatish Balay char *t_p; 280d5b485abSSatish Balay 281d5b485abSSatish Balay len = (imax+1)*(sizeof(char *) + sizeof(int *) + sizeof(int)) + 282df36731bSBarry Smith (Mbs+1)*imax*sizeof(int) + (Mbs/BITSPERBYTE+1)*imax*sizeof(char); 283d5b485abSSatish Balay table = (char **)PetscMalloc(len); CHKPTRQ(table); 284d5b485abSSatish Balay data = (int **)(table + imax); 285d5b485abSSatish Balay data[0] = (int *)(data + imax); 286df36731bSBarry Smith isz = data[0] + (Mbs+1)*imax; 287d5b485abSSatish Balay table[0] = (char *)(isz + imax); 288d5b485abSSatish Balay d_p = data[0]; t_p = table[0]; 289d5b485abSSatish Balay for (i=1; i<imax; i++) { 290df36731bSBarry Smith table[i] = t_p + (Mbs/BITSPERBYTE+1)*i; 291df36731bSBarry Smith data[i] = d_p + (Mbs+1)*i; 292d5b485abSSatish Balay } 293d5b485abSSatish Balay } 294df36731bSBarry Smith PetscMemzero(*table,(Mbs/BITSPERBYTE+1)*imax); 295d5b485abSSatish Balay PetscMemzero(isz,imax*sizeof(int)); 296d5b485abSSatish Balay 297d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 298d5b485abSSatish Balay { 299d5b485abSSatish Balay int n_i,*data_i,isz_i,*outdat_j,ctr_j; 300d5b485abSSatish Balay char *table_i; 301d5b485abSSatish Balay 302d5b485abSSatish Balay for (i=0; i<imax; i++) { 303d5b485abSSatish Balay PetscMemzero(ctr,size*sizeof(int)); 304d5b485abSSatish Balay n_i = n[i]; 305d5b485abSSatish Balay table_i = table[i]; 306d5b485abSSatish Balay idx_i = idx[i]; 307d5b485abSSatish Balay data_i = data[i]; 308d5b485abSSatish Balay isz_i = isz[i]; 309d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 310d5b485abSSatish Balay row = idx_i[j]; 311d5b485abSSatish Balay proc = rtable[row]; 312d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 313d5b485abSSatish Balay ctr[proc]++; 314d5b485abSSatish Balay *ptr[proc] = row; 315d5b485abSSatish Balay ptr[proc]++; 316d5b485abSSatish Balay } 317d5b485abSSatish Balay else { /* Update the local table */ 318d5b485abSSatish Balay if (!BT_LOOKUP(table_i,row)) { data_i[isz_i++] = row;} 319d5b485abSSatish Balay } 320d5b485abSSatish Balay } 321d5b485abSSatish Balay /* Update the headers for the current IS */ 322d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 323d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 324d5b485abSSatish Balay outdat_j = outdat[j]; 325d5b485abSSatish Balay k = ++outdat_j[0]; 326d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 327d5b485abSSatish Balay outdat_j[2*k-1] = i; 328d5b485abSSatish Balay } 329d5b485abSSatish Balay } 330d5b485abSSatish Balay isz[i] = isz_i; 331d5b485abSSatish Balay } 332d5b485abSSatish Balay } 333d5b485abSSatish Balay 334d5b485abSSatish Balay 335d5b485abSSatish Balay 336d5b485abSSatish Balay /* Now post the sends */ 337d5b485abSSatish Balay s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 338d5b485abSSatish Balay CHKPTRQ(s_waits1); 339d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 340d5b485abSSatish Balay j = pa[i]; 341d5b485abSSatish Balay MPI_Isend(outdat[j], w1[j], MPI_INT, j, tag, comm, s_waits1+i); 342d5b485abSSatish Balay } 343d5b485abSSatish Balay 344d5b485abSSatish Balay /* No longer need the original indices*/ 345d5b485abSSatish Balay for (i=0; i<imax; ++i) { 346d5b485abSSatish Balay ierr = ISRestoreIndices(is[i], idx+i); CHKERRQ(ierr); 347d5b485abSSatish Balay } 348d5b485abSSatish Balay PetscFree(idx); 349d5b485abSSatish Balay 350d5b485abSSatish Balay for (i=0; i<imax; ++i) { 351d5b485abSSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 352d5b485abSSatish Balay } 353d5b485abSSatish Balay 354d5b485abSSatish Balay /* Do Local work*/ 355df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 356d5b485abSSatish Balay 357d5b485abSSatish Balay /* Receive messages*/ 358d5b485abSSatish Balay { 359d5b485abSSatish Balay int index; 360d5b485abSSatish Balay 361d5b485abSSatish Balay recv_status = (MPI_Status *) PetscMalloc( (nrqr+1)*sizeof(MPI_Status) ); 362d5b485abSSatish Balay CHKPTRQ(recv_status); 363d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 364d5b485abSSatish Balay MPI_Waitany(nrqr, r_waits1, &index, recv_status+i); 365d5b485abSSatish Balay } 366d5b485abSSatish Balay 367d5b485abSSatish Balay s_status = (MPI_Status *) PetscMalloc( (nrqs+1)*sizeof(MPI_Status) ); 368d5b485abSSatish Balay CHKPTRQ(s_status); 369d5b485abSSatish Balay MPI_Waitall(nrqs,s_waits1,s_status); 370d5b485abSSatish Balay } 371d5b485abSSatish Balay 372d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 373d5b485abSSatish Balay PetscFree(outdat); 374d5b485abSSatish Balay PetscFree(w1); 375d5b485abSSatish Balay 376d5b485abSSatish Balay xdata = (int **)PetscMalloc((nrqr+1)*sizeof(int *)); CHKPTRQ(xdata); 377d5b485abSSatish Balay isz1 = (int *)PetscMalloc((nrqr+1)*sizeof(int)); CHKPTRQ(isz1); 378df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 379d5b485abSSatish Balay PetscFree(rbuf); 380d5b485abSSatish Balay 381d5b485abSSatish Balay /* Send the data back*/ 382d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 383d5b485abSSatish Balay { 384d5b485abSSatish Balay int *rw1, *rw2; 385d5b485abSSatish Balay 386d5b485abSSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 387d5b485abSSatish Balay PetscMemzero(rw1,2*size*sizeof(int)); 388d5b485abSSatish Balay rw2 = rw1+size; 389d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 390d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 391d5b485abSSatish Balay rw1[proc] = isz1[i]; 392d5b485abSSatish Balay } 393d5b485abSSatish Balay 394d5b485abSSatish Balay MPI_Allreduce(rw1, rw2, size, MPI_INT, MPI_MAX, comm); 395d5b485abSSatish Balay bsz1 = rw2[rank]; 396d5b485abSSatish Balay PetscFree(rw1); 397d5b485abSSatish Balay } 398d5b485abSSatish Balay 399d5b485abSSatish Balay /* Allocate buffers*/ 400d5b485abSSatish Balay 401d5b485abSSatish Balay /* Allocate memory for recv buffers. Prob none if nrqr = 0 ???? */ 402d5b485abSSatish Balay len = (nrqs+1)*sizeof(int*) + nrqs*bsz1*sizeof(int); 403d5b485abSSatish Balay rbuf2 = (int**) PetscMalloc(len); CHKPTRQ(rbuf2); 404d5b485abSSatish Balay rbuf2[0] = (int *) (rbuf2 + nrqs); 405d5b485abSSatish Balay for (i=1; i<nrqs; ++i) rbuf2[i] = rbuf2[i-1] + bsz1; 406d5b485abSSatish Balay 407d5b485abSSatish Balay /* Post the receives */ 408d5b485abSSatish Balay r_waits2 = (MPI_Request *)PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 409d5b485abSSatish Balay CHKPTRQ(r_waits2); 410d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 411d5b485abSSatish Balay MPI_Irecv(rbuf2[i], bsz1, MPI_INT, MPI_ANY_SOURCE, tag, comm, r_waits2+i); 412d5b485abSSatish Balay } 413d5b485abSSatish Balay 414d5b485abSSatish Balay /* Now post the sends */ 415d5b485abSSatish Balay s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 416d5b485abSSatish Balay CHKPTRQ(s_waits2); 417d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 418d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 419d5b485abSSatish Balay MPI_Isend( xdata[i], isz1[i], MPI_INT, j, tag, comm, s_waits2+i); 420d5b485abSSatish Balay } 421d5b485abSSatish Balay 422d5b485abSSatish Balay /* receive work done on other processors*/ 423d5b485abSSatish Balay { 424d5b485abSSatish Balay int index, is_no, ct1, max,*rbuf2_i,isz_i,*data_i,jmax; 425d5b485abSSatish Balay char *table_i; 426d5b485abSSatish Balay MPI_Status *status2; 427d5b485abSSatish Balay 428d5b485abSSatish Balay status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(status2); 429d5b485abSSatish Balay 430d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 431d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits2, &index, status2+i); 432d5b485abSSatish Balay /* Process the message*/ 433d5b485abSSatish Balay rbuf2_i = rbuf2[index]; 434d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 435d5b485abSSatish Balay jmax = rbuf2[index][0]; 436d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 437d5b485abSSatish Balay max = rbuf2_i[2*j]; 438d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 439d5b485abSSatish Balay isz_i = isz[is_no]; 440d5b485abSSatish Balay data_i = data[is_no]; 441d5b485abSSatish Balay table_i = table[is_no]; 442d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 443d5b485abSSatish Balay row = rbuf2_i[ct1]; 444d5b485abSSatish Balay if (!BT_LOOKUP(table_i,row)) { data_i[isz_i++] = row;} 445d5b485abSSatish Balay } 446d5b485abSSatish Balay isz[is_no] = isz_i; 447d5b485abSSatish Balay } 448d5b485abSSatish Balay } 449d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits2,status2); 450d5b485abSSatish Balay PetscFree(status2); 451d5b485abSSatish Balay } 452d5b485abSSatish Balay 453d5b485abSSatish Balay for (i=0; i<imax; ++i) { 454d5b485abSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz[i], data[i], is+i); CHKERRQ(ierr); 455d5b485abSSatish Balay } 456d5b485abSSatish Balay 457d5b485abSSatish Balay PetscFree(pa); 458d5b485abSSatish Balay PetscFree(rbuf2); 459d5b485abSSatish Balay PetscFree(s_waits1); 460d5b485abSSatish Balay PetscFree(r_waits1); 461d5b485abSSatish Balay PetscFree(s_waits2); 462d5b485abSSatish Balay PetscFree(r_waits2); 463d5b485abSSatish Balay PetscFree(table); 464d5b485abSSatish Balay PetscFree(s_status); 465d5b485abSSatish Balay PetscFree(recv_status); 466d5b485abSSatish Balay PetscFree(xdata[0]); 467d5b485abSSatish Balay PetscFree(xdata); 468d5b485abSSatish Balay PetscFree(isz1); 469d5b485abSSatish Balay return 0; 470d5b485abSSatish Balay } 471d5b485abSSatish Balay 472d5b485abSSatish Balay /* 473df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 474d5b485abSSatish Balay the work on the local processor. 475d5b485abSSatish Balay 476d5b485abSSatish Balay Inputs: 477df36731bSBarry Smith C - MAT_MPIBAIJ; 478d5b485abSSatish Balay imax - total no of index sets processed at a time; 479df36731bSBarry Smith table - an array of char - size = Mbs bits. 480d5b485abSSatish Balay 481d5b485abSSatish Balay Output: 482d5b485abSSatish Balay isz - array containing the count of the solution elements correspondign 483d5b485abSSatish Balay to each index set; 484d5b485abSSatish Balay data - pointer to the solutions 485d5b485abSSatish Balay */ 486df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,char **table,int *isz, 487d5b485abSSatish Balay int **data) 488d5b485abSSatish Balay { 489df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 490d5b485abSSatish Balay Mat A = c->A, B = c->B; 491df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 492df36731bSBarry Smith int start, end, val, max, rstart,cstart,*ai, *aj; 493d5b485abSSatish Balay int *bi, *bj, *garray, i, j, k, row,*data_i,isz_i; 494d5b485abSSatish Balay char *table_i; 495d5b485abSSatish Balay 496d5b485abSSatish Balay rstart = c->rstart; 497d5b485abSSatish Balay cstart = c->cstart; 498d5b485abSSatish Balay ai = a->i; 499df36731bSBarry Smith aj = a->j; 500d5b485abSSatish Balay bi = b->i; 501df36731bSBarry Smith bj = b->j; 502d5b485abSSatish Balay garray = c->garray; 503d5b485abSSatish Balay 504d5b485abSSatish Balay 505d5b485abSSatish Balay for (i=0; i<imax; i++) { 506d5b485abSSatish Balay data_i = data[i]; 507d5b485abSSatish Balay table_i = table[i]; 508d5b485abSSatish Balay isz_i = isz[i]; 509d5b485abSSatish Balay for (j=0, max=isz[i]; j<max; j++) { 510d5b485abSSatish Balay row = data_i[j] - rstart; 511d5b485abSSatish Balay start = ai[row]; 512d5b485abSSatish Balay end = ai[row+1]; 513d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 514df36731bSBarry Smith val = aj[k] + cstart; 515d5b485abSSatish Balay if (!BT_LOOKUP(table_i,val)) { data_i[isz_i++] = val;} 516d5b485abSSatish Balay } 517d5b485abSSatish Balay start = bi[row]; 518d5b485abSSatish Balay end = bi[row+1]; 519d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 520df36731bSBarry Smith val = garray[bj[k]]; 521d5b485abSSatish Balay if (!BT_LOOKUP(table_i,val)) { data_i[isz_i++] = val;} 522d5b485abSSatish Balay } 523d5b485abSSatish Balay } 524d5b485abSSatish Balay isz[i] = isz_i; 525d5b485abSSatish Balay } 526d5b485abSSatish Balay return 0; 527d5b485abSSatish Balay } 528d5b485abSSatish Balay /* 529df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 530d5b485abSSatish Balay and return the output 531d5b485abSSatish Balay 532d5b485abSSatish Balay Input: 533d5b485abSSatish Balay C - the matrix 534d5b485abSSatish Balay nrqr - no of messages being processed. 535d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 536d5b485abSSatish Balay 537d5b485abSSatish Balay Output: 538d5b485abSSatish Balay xdata - array of messages to be sent back 539d5b485abSSatish Balay isz1 - size of each message 540d5b485abSSatish Balay 541d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 542d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 543d5b485abSSatish Balay rather then all previous rows as it is now where a single large chunck of 544d5b485abSSatish Balay memory is used. 545d5b485abSSatish Balay 546d5b485abSSatish Balay */ 547df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf, 548d5b485abSSatish Balay int **xdata, int * isz1) 549d5b485abSSatish Balay { 550df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 551d5b485abSSatish Balay Mat A = c->A, B = c->B; 552df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 553df36731bSBarry Smith int rstart,cstart,*ai, *aj, *bi, *bj, *garray, i, j, k; 554d5b485abSSatish Balay int row,total_sz,ct, ct1, ct2, ct3,mem_estimate, oct2, l, start, end; 555df36731bSBarry Smith int val, max1, max2, rank, Mbs, no_malloc =0, *tmp, new_estimate, ctr; 556d5b485abSSatish Balay int *rbuf_i,kmax,rbuf_0; 557d5b485abSSatish Balay char *xtable; 558d5b485abSSatish Balay 559d5b485abSSatish Balay rank = c->rank; 560df36731bSBarry Smith Mbs = c->Mbs; 561d5b485abSSatish Balay rstart = c->rstart; 562d5b485abSSatish Balay cstart = c->cstart; 563d5b485abSSatish Balay ai = a->i; 564df36731bSBarry Smith aj = a->j; 565d5b485abSSatish Balay bi = b->i; 566df36731bSBarry Smith bj = b->j; 567d5b485abSSatish Balay garray = c->garray; 568d5b485abSSatish Balay 569d5b485abSSatish Balay 570d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 571d5b485abSSatish Balay rbuf_i = rbuf[i]; 572d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 573d5b485abSSatish Balay ct += rbuf_0; 574d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 575d5b485abSSatish Balay } 576d5b485abSSatish Balay 577df36731bSBarry Smith max1 = ct*(a->nz +b->nz)/c->Mbs; 578d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 579d5b485abSSatish Balay xdata[0] = (int *)PetscMalloc(mem_estimate*sizeof(int)); CHKPTRQ(xdata[0]); 580d5b485abSSatish Balay ++no_malloc; 581df36731bSBarry Smith xtable = (char *)PetscMalloc((Mbs/BITSPERBYTE+1)*sizeof(char)); CHKPTRQ(xtable); 582d5b485abSSatish Balay PetscMemzero(isz1,nrqr*sizeof(int)); 583d5b485abSSatish Balay 584d5b485abSSatish Balay ct3 = 0; 585d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 586d5b485abSSatish Balay rbuf_i = rbuf[i]; 587d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 588d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 589d5b485abSSatish Balay ct2 = ct1; 590d5b485abSSatish Balay ct3 += ct1; 591d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 592df36731bSBarry Smith PetscMemzero(xtable,(Mbs/BITSPERBYTE+1)*sizeof(char)); 593d5b485abSSatish Balay oct2 = ct2; 594d5b485abSSatish Balay kmax = rbuf_i[2*j]; 595d5b485abSSatish Balay for (k=0; k<kmax; k++, ct1++) { 596d5b485abSSatish Balay row = rbuf_i[ct1]; 597d5b485abSSatish Balay if (!BT_LOOKUP(xtable,row)) { 598d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 599d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 600d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 601d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 602d5b485abSSatish Balay PetscFree(xdata[0]); 603d5b485abSSatish Balay xdata[0] = tmp; 604d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 605d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 606d5b485abSSatish Balay } 607d5b485abSSatish Balay xdata[i][ct2++] = row; 608d5b485abSSatish Balay ct3++; 609d5b485abSSatish Balay } 610d5b485abSSatish Balay } 611d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 612d5b485abSSatish Balay row = xdata[i][k] - rstart; 613d5b485abSSatish Balay start = ai[row]; 614d5b485abSSatish Balay end = ai[row+1]; 615d5b485abSSatish Balay for (l=start; l<end; l++) { 616df36731bSBarry Smith val = aj[l] + cstart; 617d5b485abSSatish Balay if (!BT_LOOKUP(xtable,val)) { 618d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 619d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 620d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 621d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 622d5b485abSSatish Balay PetscFree(xdata[0]); 623d5b485abSSatish Balay xdata[0] = tmp; 624d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 625d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 626d5b485abSSatish Balay } 627d5b485abSSatish Balay xdata[i][ct2++] = val; 628d5b485abSSatish Balay ct3++; 629d5b485abSSatish Balay } 630d5b485abSSatish Balay } 631d5b485abSSatish Balay start = bi[row]; 632d5b485abSSatish Balay end = bi[row+1]; 633d5b485abSSatish Balay for (l=start; l<end; l++) { 634df36731bSBarry Smith val = garray[bj[l]]; 635d5b485abSSatish Balay if (!BT_LOOKUP(xtable,val)) { 636d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 637d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 638d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 639d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 640d5b485abSSatish Balay PetscFree(xdata[0]); 641d5b485abSSatish Balay xdata[0] = tmp; 642d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 643d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 644d5b485abSSatish Balay } 645d5b485abSSatish Balay xdata[i][ct2++] = val; 646d5b485abSSatish Balay ct3++; 647d5b485abSSatish Balay } 648d5b485abSSatish Balay } 649d5b485abSSatish Balay } 650d5b485abSSatish Balay /* Update the header*/ 651d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 652d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 653d5b485abSSatish Balay } 654d5b485abSSatish Balay xdata[i][0] = rbuf_0; 655d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 656d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 657d5b485abSSatish Balay } 658d5b485abSSatish Balay PetscFree(xtable); 659df36731bSBarry Smith PLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d bytes, no of mallocs = %d\n",rank,mem_estimate, ct3,no_malloc); 660d5b485abSSatish Balay return 0; 661d5b485abSSatish Balay } 662d5b485abSSatish Balay 663a2feddc5SSatish Balay static int MatGetSubMatrices_MPIBAIJ_local(Mat,int,IS *,IS *,MatGetSubMatrixCall,Mat **); 664a2feddc5SSatish Balay 665df36731bSBarry Smith int MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,IS *isrow,IS *iscol, 666d5b485abSSatish Balay MatGetSubMatrixCall scall,Mat **submat) 667d5b485abSSatish Balay { 66836f4e84dSSatish Balay int i,ierr; 66936f4e84dSSatish Balay IS *isrow_new,*iscol_new; 670a2feddc5SSatish Balay 671a2feddc5SSatish Balay /* The compression and expansion should be avoided. Dos'nt point 672a2feddc5SSatish Balay out errors might change the indices hence buggey */ 673a2feddc5SSatish Balay 67436f4e84dSSatish Balay isrow_new = (IS *)PetscMalloc(2*ismax*sizeof(IS)); CHKPTRQ(isrow_new); 67536f4e84dSSatish Balay iscol_new = isrow_new + ismax; 676e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C, ismax, isrow,isrow_new); CHKERRQ(ierr); 677e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C, ismax, iscol,iscol_new); CHKERRQ(ierr); 67836f4e84dSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,ismax,isrow_new,iscol_new,scall,submat); CHKERRQ(ierr); 67936f4e84dSSatish Balay 68036f4e84dSSatish Balay for (i=0; i<ismax; i++) { 68136f4e84dSSatish Balay ISDestroy(isrow_new[i]); 68236f4e84dSSatish Balay ISDestroy(iscol_new[i]); 68336f4e84dSSatish Balay } 68436f4e84dSSatish Balay PetscFree(isrow_new); 685a2feddc5SSatish Balay return 0; 686a2feddc5SSatish Balay } 687a2feddc5SSatish Balay 688a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 689a2feddc5SSatish Balay static int MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,IS *isrow,IS *iscol, 690a2feddc5SSatish Balay MatGetSubMatrixCall scall,Mat **submat) 691a2feddc5SSatish Balay { 692df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 693d5b485abSSatish Balay Mat A = c->A,*submats = *submat; 694df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data, *b = (Mat_SeqBAIJ*)c->B->data, *mat; 695d5b485abSSatish Balay int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size; 696a2feddc5SSatish Balay int **sbuf1,**sbuf2, rank, Mbs,i,j,k,l,ct1,ct2,ierr, **rbuf1,row,proc; 69734e6ae18SSatish Balay int nrqs, msz, **ptr,index,*req_size,*ctr,*pa,*tmp,tcol,bsz,nrqr; 698df36731bSBarry Smith int **rbuf3,*req_source,**sbuf_aj, **rbuf2, max1,max2,**rmap; 699d5b485abSSatish Balay int **cmap,**lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 700d5b485abSSatish Balay int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*cmap_i,*lens_i; 70136f4e84dSSatish Balay int *rmap_i,bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA, *cworkB; 70236f4e84dSSatish Balay int cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 70334e6ae18SSatish Balay int *bmap = c->garray,ctmp,rstart=c->rstart,tag0,tag1,tag2,tag3; 704d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 705d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 706d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 707d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 708d5b485abSSatish Balay MPI_Comm comm; 70936f4e84dSSatish Balay Scalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 71036f4e84dSSatish Balay Scalar *a_a=a->a,*b_a=b->a; 711d5b485abSSatish Balay 712d5b485abSSatish Balay comm = C->comm; 71334e6ae18SSatish Balay tag0 = C->tag; 714d5b485abSSatish Balay size = c->size; 715d5b485abSSatish Balay rank = c->rank; 716a2feddc5SSatish Balay Mbs = c->Mbs; 717d5b485abSSatish Balay 71834e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 71934e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1); CHKERRQ(ierr); 72034e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2); CHKERRQ(ierr); 72134e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3); CHKERRQ(ierr); 72234e6ae18SSatish Balay 723d5b485abSSatish Balay /* Check if the col indices are sorted */ 724d5b485abSSatish Balay for (i=0; i<ismax; i++) { 725d5b485abSSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j); 726df36731bSBarry Smith if (!j) SETERRQ(1,"MatGetSubmatrices_MPIBAIJ:IS is not sorted"); 727d5b485abSSatish Balay } 728d5b485abSSatish Balay 729a2feddc5SSatish Balay len = (2*ismax+1)*(sizeof(int *) + sizeof(int)) + (Mbs+1)*sizeof(int); 730d5b485abSSatish Balay irow = (int **)PetscMalloc(len); CHKPTRQ(irow); 731d5b485abSSatish Balay icol = irow + ismax; 732d5b485abSSatish Balay nrow = (int *) (icol + ismax); 733d5b485abSSatish Balay ncol = nrow + ismax; 734d5b485abSSatish Balay rtable = ncol + ismax; 735d5b485abSSatish Balay 736d5b485abSSatish Balay for (i=0; i<ismax; i++) { 737d5b485abSSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]); CHKERRQ(ierr); 738d5b485abSSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]); CHKERRQ(ierr); 739d5b485abSSatish Balay ierr = ISGetSize(isrow[i],&nrow[i]); CHKERRQ(ierr); 740d5b485abSSatish Balay ierr = ISGetSize(iscol[i],&ncol[i]); CHKERRQ(ierr); 741d5b485abSSatish Balay } 742d5b485abSSatish Balay 743d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 744d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 745d5b485abSSatish Balay jmax = c->rowners[i+1]; 746d5b485abSSatish Balay for (; j<jmax; j++) { 747d5b485abSSatish Balay rtable[j] = i; 748d5b485abSSatish Balay } 749d5b485abSSatish Balay } 750d5b485abSSatish Balay 751d5b485abSSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 752d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 753d5b485abSSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */ 754d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 755d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 756d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 757d5b485abSSatish Balay PetscMemzero(w1,size*3*sizeof(int)); /* initialise work vector*/ 758d5b485abSSatish Balay for (i=0; i<ismax; i++) { 759d5b485abSSatish Balay PetscMemzero(w4,size*sizeof(int)); /* initialise work vector*/ 760d5b485abSSatish Balay jmax = nrow[i]; 761d5b485abSSatish Balay irow_i = irow[i]; 762d5b485abSSatish Balay for (j=0; j<jmax; j++) { 763d5b485abSSatish Balay row = irow_i[j]; 764d5b485abSSatish Balay proc = rtable[row]; 765d5b485abSSatish Balay w4[proc]++; 766d5b485abSSatish Balay } 767d5b485abSSatish Balay for (j=0; j<size; j++) { 768d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 769d5b485abSSatish Balay } 770d5b485abSSatish Balay } 771d5b485abSSatish Balay 772d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 773e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 774d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 775d5b485abSSatish Balay w3[rank] = 0; 776d5b485abSSatish Balay for (i=0; i<size; i++) { 777d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 778d5b485abSSatish Balay } 779d5b485abSSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/ 780d5b485abSSatish Balay for (i=0, j=0; i<size; i++) { 781d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 782d5b485abSSatish Balay } 783d5b485abSSatish Balay 784d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 785d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 786d5b485abSSatish Balay j = pa[i]; 787d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 788d5b485abSSatish Balay msz += w1[j]; 789d5b485abSSatish Balay } 790d5b485abSSatish Balay /* Do a global reduction to determine how many messages to expect*/ 791d5b485abSSatish Balay { 792d5b485abSSatish Balay int *rw1, *rw2; 793d5b485abSSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 794d5b485abSSatish Balay rw2 = rw1+size; 795d5b485abSSatish Balay MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm); 796d5b485abSSatish Balay bsz = rw1[rank]; 797d5b485abSSatish Balay MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm); 798d5b485abSSatish Balay nrqr = rw2[rank]; 799d5b485abSSatish Balay PetscFree(rw1); 800d5b485abSSatish Balay } 801d5b485abSSatish Balay 802d5b485abSSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 803d5b485abSSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 804d5b485abSSatish Balay rbuf1 = (int**) PetscMalloc(len); CHKPTRQ(rbuf1); 805d5b485abSSatish Balay rbuf1[0] = (int *) (rbuf1 + nrqr); 806d5b485abSSatish Balay for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz; 807d5b485abSSatish Balay 808d5b485abSSatish Balay /* Post the receives */ 809d5b485abSSatish Balay r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 810d5b485abSSatish Balay CHKPTRQ(r_waits1); 811d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 81234e6ae18SSatish Balay MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i); 813d5b485abSSatish Balay } 814d5b485abSSatish Balay 815d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 816d5b485abSSatish Balay len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 817d5b485abSSatish Balay sbuf1 = (int **)PetscMalloc(len); CHKPTRQ(sbuf1); 818d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 819d5b485abSSatish Balay PetscMemzero(sbuf1,2*size*sizeof(int*)); 820d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 821d5b485abSSatish Balay tmp = (int *) (ptr + size); 822d5b485abSSatish Balay ctr = tmp + 2*msz; 823d5b485abSSatish Balay 824d5b485abSSatish Balay { 82536f4e84dSSatish Balay 826d5b485abSSatish Balay int *iptr = tmp,ict = 0; 827d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 828d5b485abSSatish Balay j = pa[i]; 829d5b485abSSatish Balay iptr += ict; 830d5b485abSSatish Balay sbuf1[j] = iptr; 831d5b485abSSatish Balay ict = w1[j]; 832d5b485abSSatish Balay } 833d5b485abSSatish Balay } 834d5b485abSSatish Balay 835d5b485abSSatish Balay /* Form the outgoing messages */ 836d5b485abSSatish Balay /* Initialise the header space */ 837d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 838d5b485abSSatish Balay j = pa[i]; 839d5b485abSSatish Balay sbuf1[j][0] = 0; 840d5b485abSSatish Balay PetscMemzero(sbuf1[j]+1, 2*w3[j]*sizeof(int)); 841d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 842d5b485abSSatish Balay } 843d5b485abSSatish Balay 844d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 845d5b485abSSatish Balay for (i=0; i<ismax; i++) { 846d5b485abSSatish Balay PetscMemzero(ctr,size*sizeof(int)); 847d5b485abSSatish Balay irow_i = irow[i]; 848d5b485abSSatish Balay jmax = nrow[i]; 849d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 850d5b485abSSatish Balay row = irow_i[j]; 851d5b485abSSatish Balay proc = rtable[row]; 852d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 853d5b485abSSatish Balay ctr[proc]++; 854d5b485abSSatish Balay *ptr[proc] = row; 855d5b485abSSatish Balay ptr[proc]++; 856d5b485abSSatish Balay } 857d5b485abSSatish Balay } 858d5b485abSSatish Balay /* Update the headers for the current IS */ 859d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 860d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 861d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 862d5b485abSSatish Balay k = ++sbuf1_j[0]; 863d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 864d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 865d5b485abSSatish Balay } 866d5b485abSSatish Balay } 867d5b485abSSatish Balay } 868d5b485abSSatish Balay 869d5b485abSSatish Balay /* Now post the sends */ 870d5b485abSSatish Balay s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 871d5b485abSSatish Balay CHKPTRQ(s_waits1); 872d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 873d5b485abSSatish Balay j = pa[i]; 874d5b485abSSatish Balay /* printf("[%d] Send Req to %d: size %d \n", rank,j, w1[j]); */ 87534e6ae18SSatish Balay MPI_Isend( sbuf1[j], w1[j], MPI_INT, j, tag0, comm, s_waits1+i); 876d5b485abSSatish Balay } 877d5b485abSSatish Balay 878d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 879d5b485abSSatish Balay r_waits2 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 880d5b485abSSatish Balay CHKPTRQ(r_waits2); 881d5b485abSSatish Balay rbuf2 = (int**)PetscMalloc((nrqs+1)*sizeof(int *));CHKPTRQ(rbuf2); 882d5b485abSSatish Balay rbuf2[0] = tmp + msz; 883d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 884d5b485abSSatish Balay j = pa[i]; 885d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 886d5b485abSSatish Balay } 887d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 888d5b485abSSatish Balay j = pa[i]; 88934e6ae18SSatish Balay MPI_Irecv( rbuf2[i], w1[j], MPI_INT, j, tag1, comm, r_waits2+i); 890d5b485abSSatish Balay } 891d5b485abSSatish Balay 892d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 893d5b485abSSatish Balay 894d5b485abSSatish Balay 895d5b485abSSatish Balay /* Receive messages*/ 896d5b485abSSatish Balay s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 897d5b485abSSatish Balay CHKPTRQ(s_waits2); 898d5b485abSSatish Balay r_status1 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 899d5b485abSSatish Balay CHKPTRQ(r_status1); 900d5b485abSSatish Balay len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*); 901d5b485abSSatish Balay sbuf2 = (int**) PetscMalloc(len);CHKPTRQ(sbuf2); 902d5b485abSSatish Balay req_size = (int *) (sbuf2 + nrqr); 903d5b485abSSatish Balay req_source = req_size + nrqr; 904d5b485abSSatish Balay 905d5b485abSSatish Balay { 906df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*) c->A->data, *sB = (Mat_SeqBAIJ*) c->B->data; 90736f4e84dSSatish Balay int *sAi = sA->i, *sBi = sB->i, id; 908d5b485abSSatish Balay int *sbuf2_i; 909d5b485abSSatish Balay 910d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 911d5b485abSSatish Balay MPI_Waitany(nrqr, r_waits1, &index, r_status1+i); 912d5b485abSSatish Balay req_size[index] = 0; 913d5b485abSSatish Balay rbuf1_i = rbuf1[index]; 914d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 915d5b485abSSatish Balay MPI_Get_count(r_status1+i,MPI_INT, &end); 916d5b485abSSatish Balay sbuf2[index] = (int *)PetscMalloc(end*sizeof(int));CHKPTRQ(sbuf2[index]); 917d5b485abSSatish Balay sbuf2_i = sbuf2[index]; 918d5b485abSSatish Balay for (j=start; j<end; j++) { 919d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 920d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 921d5b485abSSatish Balay sbuf2_i[j] = ncols; 922d5b485abSSatish Balay req_size[index] += ncols; 923d5b485abSSatish Balay } 924d5b485abSSatish Balay req_source[index] = r_status1[i].MPI_SOURCE; 925d5b485abSSatish Balay /* form the header */ 926d5b485abSSatish Balay sbuf2_i[0] = req_size[index]; 927d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 92834e6ae18SSatish Balay MPI_Isend(sbuf2_i,end,MPI_INT,req_source[index],tag1,comm,s_waits2+i); 929d5b485abSSatish Balay } 930d5b485abSSatish Balay } 931d5b485abSSatish Balay PetscFree(r_status1); PetscFree(r_waits1); 932d5b485abSSatish Balay 933d5b485abSSatish Balay /* recv buffer sizes */ 934d5b485abSSatish Balay /* Receive messages*/ 935d5b485abSSatish Balay 936d5b485abSSatish Balay rbuf3 = (int**)PetscMalloc((nrqs+1)*sizeof(int*)); CHKPTRQ(rbuf3); 937d5b485abSSatish Balay rbuf4 = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar*));CHKPTRQ(rbuf4); 938d5b485abSSatish Balay r_waits3 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 939d5b485abSSatish Balay CHKPTRQ(r_waits3); 940d5b485abSSatish Balay r_waits4 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 941d5b485abSSatish Balay CHKPTRQ(r_waits4); 942d5b485abSSatish Balay r_status2 = (MPI_Status *) PetscMalloc( (nrqs+1)*sizeof(MPI_Status) ); 943d5b485abSSatish Balay CHKPTRQ(r_status2); 944d5b485abSSatish Balay 945d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 946d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits2, &index, r_status2+i); 947d5b485abSSatish Balay rbuf3[index] = (int *)PetscMalloc(rbuf2[index][0]*sizeof(int)); 948d5b485abSSatish Balay CHKPTRQ(rbuf3[index]); 949a2feddc5SSatish Balay rbuf4[index] = (Scalar *)PetscMalloc(rbuf2[index][0]*bs2*sizeof(Scalar)); 950d5b485abSSatish Balay CHKPTRQ(rbuf4[index]); 951d5b485abSSatish Balay MPI_Irecv(rbuf3[index],rbuf2[index][0], MPI_INT, 95234e6ae18SSatish Balay r_status2[i].MPI_SOURCE, tag2, comm, r_waits3+index); 953a2feddc5SSatish Balay MPI_Irecv(rbuf4[index],rbuf2[index][0]*bs2, MPIU_SCALAR, 95434e6ae18SSatish Balay r_status2[i].MPI_SOURCE, tag3, comm, r_waits4+index); 955d5b485abSSatish Balay } 956d5b485abSSatish Balay PetscFree(r_status2); PetscFree(r_waits2); 957d5b485abSSatish Balay 958d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 959d5b485abSSatish Balay s_status1 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status)); 960d5b485abSSatish Balay CHKPTRQ(s_status1); 961d5b485abSSatish Balay s_status2 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 962d5b485abSSatish Balay CHKPTRQ(s_status2); 963d5b485abSSatish Balay 964d5b485abSSatish Balay MPI_Waitall(nrqs,s_waits1,s_status1); 965d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits2,s_status2); 966d5b485abSSatish Balay PetscFree(s_status1); PetscFree(s_status2); 967d5b485abSSatish Balay PetscFree(s_waits1); PetscFree(s_waits2); 968d5b485abSSatish Balay 969d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 970d5b485abSSatish Balay sbuf_aj = (int **)PetscMalloc((nrqr+1)*sizeof(int *));CHKPTRQ(sbuf_aj); 971d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 972d5b485abSSatish Balay sbuf_aj[0] = (int*) PetscMalloc((j+1)*sizeof(int)); CHKPTRQ(sbuf_aj[0]); 973d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 974d5b485abSSatish Balay 975d5b485abSSatish Balay s_waits3 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 976d5b485abSSatish Balay CHKPTRQ(s_waits3); 977d5b485abSSatish Balay { 978d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 979d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 980d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 981d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 982d5b485abSSatish Balay ct2 = 0; 983d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 984d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 985d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 986e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 987d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 988d5b485abSSatish Balay ncols = nzA + nzB; 989d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 990d5b485abSSatish Balay 991d5b485abSSatish Balay /* load the column indices for this row into cols*/ 992d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 993d5b485abSSatish Balay for (l=0; l<nzB; l++) { 994d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 995d5b485abSSatish Balay else break; 996d5b485abSSatish Balay } 997d5b485abSSatish Balay imark = l; 998d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 999d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1000d5b485abSSatish Balay ct2 += ncols; 1001d5b485abSSatish Balay } 1002d5b485abSSatish Balay } 100334e6ae18SSatish Balay MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i); 1004d5b485abSSatish Balay } 1005d5b485abSSatish Balay } 1006d5b485abSSatish Balay r_status3 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status)); 1007d5b485abSSatish Balay CHKPTRQ(r_status3); 1008d5b485abSSatish Balay s_status3 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 1009d5b485abSSatish Balay CHKPTRQ(s_status3); 1010d5b485abSSatish Balay 1011d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 1012d5b485abSSatish Balay sbuf_aa = (Scalar **)PetscMalloc((nrqr+1)*sizeof(Scalar *));CHKPTRQ(sbuf_aa); 1013d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1014a2feddc5SSatish Balay sbuf_aa[0] = (Scalar*) PetscMalloc((j+1)*bs2*sizeof(Scalar));CHKPTRQ(sbuf_aa[0]); 1015a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1016d5b485abSSatish Balay 1017d5b485abSSatish Balay s_waits4 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 1018d5b485abSSatish Balay CHKPTRQ(s_waits4); 1019d5b485abSSatish Balay { 1020d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1021d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1022d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 1023d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 1024d5b485abSSatish Balay ct2 = 0; 1025d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1026d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 1027d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1028e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1029d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1030d5b485abSSatish Balay ncols = nzA + nzB; 1031d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1032a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 1033d5b485abSSatish Balay 1034d5b485abSSatish Balay /* load the column values for this row into vals*/ 10355b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 1036d5b485abSSatish Balay for (l=0; l<nzB; l++) { 1037a2feddc5SSatish Balay if ((bmap[cworkB[l]]) < cstart) 1038a2feddc5SSatish Balay PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(Scalar)); 1039d5b485abSSatish Balay else break; 1040d5b485abSSatish Balay } 1041d5b485abSSatish Balay imark = l; 1042a2feddc5SSatish Balay for (l=0; l<nzA; l++) 1043a2feddc5SSatish Balay PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(Scalar)); 1044a2feddc5SSatish Balay for (l=imark; l<nzB; l++) 1045a2feddc5SSatish Balay PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(Scalar)); 1046d5b485abSSatish Balay ct2 += ncols; 1047d5b485abSSatish Balay } 1048d5b485abSSatish Balay } 104934e6ae18SSatish Balay MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i); 1050d5b485abSSatish Balay } 1051d5b485abSSatish Balay } 1052d5b485abSSatish Balay r_status4 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status)); 1053d5b485abSSatish Balay CHKPTRQ(r_status4); 1054d5b485abSSatish Balay s_status4 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 1055d5b485abSSatish Balay CHKPTRQ(s_status4); 1056d5b485abSSatish Balay PetscFree(rbuf1); 1057d5b485abSSatish Balay 1058d5b485abSSatish Balay /* Form the matrix */ 1059d5b485abSSatish Balay /* create col map */ 1060d5b485abSSatish Balay { 1061d5b485abSSatish Balay int *icol_i; 1062d5b485abSSatish Balay 1063a2feddc5SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->Nbs*sizeof(int); 1064d5b485abSSatish Balay cmap = (int **)PetscMalloc(len); CHKPTRQ(cmap); 1065d5b485abSSatish Balay cmap[0] = (int *)(cmap + ismax); 1066a2feddc5SSatish Balay PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int)); 1067a2feddc5SSatish Balay for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; } 1068d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1069d5b485abSSatish Balay jmax = ncol[i]; 1070d5b485abSSatish Balay icol_i = icol[i]; 1071d5b485abSSatish Balay cmap_i = cmap[i]; 1072d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1073d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1074d5b485abSSatish Balay } 1075d5b485abSSatish Balay } 1076d5b485abSSatish Balay } 1077d5b485abSSatish Balay 1078d5b485abSSatish Balay 1079d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1080d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1081d5b485abSSatish Balay len = (1+ismax)*sizeof(int *) + j*sizeof(int); 1082d5b485abSSatish Balay lens = (int **)PetscMalloc(len); CHKPTRQ(lens); 1083d5b485abSSatish Balay lens[0] = (int *)(lens + ismax); 1084d5b485abSSatish Balay PetscMemzero(lens[0], j*sizeof(int)); 1085d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1086d5b485abSSatish Balay 1087d5b485abSSatish Balay /* Update lens from local data */ 1088d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1089d5b485abSSatish Balay jmax = nrow[i]; 1090d5b485abSSatish Balay cmap_i = cmap[i]; 1091d5b485abSSatish Balay irow_i = irow[i]; 1092d5b485abSSatish Balay lens_i = lens[i]; 1093d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1094d5b485abSSatish Balay row = irow_i[j]; 1095d5b485abSSatish Balay proc = rtable[row]; 1096d5b485abSSatish Balay if (proc == rank) { 109736f4e84dSSatish Balay /* Get indices from matA and then from matB */ 109836f4e84dSSatish Balay row = row - rstart; 109936f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 110036f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 110136f4e84dSSatish Balay for (k=0; k<nzA; k++) 110236f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++;} 110336f4e84dSSatish Balay for (k=0; k<nzB; k++) 110436f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++;} 1105d5b485abSSatish Balay } 1106d5b485abSSatish Balay } 1107a2feddc5SSatish Balay } 1108d5b485abSSatish Balay 1109d5b485abSSatish Balay /* Create row map*/ 1110a2feddc5SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->Mbs*sizeof(int); 1111d5b485abSSatish Balay rmap = (int **)PetscMalloc(len); CHKPTRQ(rmap); 1112d5b485abSSatish Balay rmap[0] = (int *)(rmap + ismax); 1113a2feddc5SSatish Balay PetscMemzero(rmap[0],ismax*c->Mbs*sizeof(int)); 1114a2feddc5SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + c->Mbs;} 1115d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1116d5b485abSSatish Balay rmap_i = rmap[i]; 1117d5b485abSSatish Balay irow_i = irow[i]; 1118d5b485abSSatish Balay jmax = nrow[i]; 1119d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1120d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1121d5b485abSSatish Balay } 1122d5b485abSSatish Balay } 1123d5b485abSSatish Balay 1124d5b485abSSatish Balay /* Update lens from offproc data */ 1125d5b485abSSatish Balay { 1126d5b485abSSatish Balay int *rbuf2_i, *rbuf3_i, *sbuf1_i; 1127d5b485abSSatish Balay 1128d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1129d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits3, &i, r_status3+tmp2); 1130d5b485abSSatish Balay index = pa[i]; 1131d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1132d5b485abSSatish Balay jmax = sbuf1_i[0]; 1133d5b485abSSatish Balay ct1 = 2*jmax+1; 1134d5b485abSSatish Balay ct2 = 0; 1135d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1136d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1137d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1138d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1139d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1140d5b485abSSatish Balay lens_i = lens[is_no]; 1141d5b485abSSatish Balay cmap_i = cmap[is_no]; 1142d5b485abSSatish Balay rmap_i = rmap[is_no]; 1143d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1144d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1145d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1146d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1147d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1148d5b485abSSatish Balay lens_i[row]++; 1149d5b485abSSatish Balay } 1150d5b485abSSatish Balay } 1151d5b485abSSatish Balay } 1152d5b485abSSatish Balay } 1153d5b485abSSatish Balay } 1154d5b485abSSatish Balay } 1155d5b485abSSatish Balay PetscFree(r_status3); PetscFree(r_waits3); 1156d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits3,s_status3); 1157d5b485abSSatish Balay PetscFree(s_status3); PetscFree(s_waits3); 1158d5b485abSSatish Balay 1159d5b485abSSatish Balay /* Create the submatrices */ 1160d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1161d5b485abSSatish Balay /* 1162d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1163d5b485abSSatish Balay */ 1164d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1165df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 116636f4e84dSSatish Balay if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) { 1167df36731bSBarry Smith SETERRQ(1,"MatGetSubmatrices_MPIBAIJ:Cannot reuse matrix. wrong size"); 1168d5b485abSSatish Balay } 116936f4e84dSSatish Balay if (PetscMemcmp(mat->ilen,lens[i], mat->mbs *sizeof(int))) { 1170df36731bSBarry Smith SETERRQ(1,"MatGetSubmatrices_MPIBAIJ:Cannot reuse matrix. wrong no of nonzeros"); 1171d5b485abSSatish Balay } 1172d5b485abSSatish Balay /* Initial matrix as if empty */ 117336f4e84dSSatish Balay PetscMemzero(mat->ilen,mat->mbs*sizeof(int)); 1174d5b485abSSatish Balay } 1175d5b485abSSatish Balay } 1176d5b485abSSatish Balay else { 1177d5b485abSSatish Balay *submat = submats = (Mat *)PetscMalloc(ismax*sizeof(Mat)); CHKPTRQ(submats); 1178d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1179a2feddc5SSatish Balay ierr = MatCreateSeqBAIJ(comm,a->bs,nrow[i]*bs,ncol[i]*bs,0,lens[i],submats+i);CHKERRQ(ierr); 1180d5b485abSSatish Balay } 1181d5b485abSSatish Balay } 1182d5b485abSSatish Balay 1183d5b485abSSatish Balay /* Assemble the matrices */ 1184d5b485abSSatish Balay /* First assemble the local rows */ 1185d5b485abSSatish Balay { 118636f4e84dSSatish Balay int ilen_row,*imat_ilen, *imat_j, *imat_i; 1187d5b485abSSatish Balay Scalar *imat_a; 1188d5b485abSSatish Balay 1189d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1190df36731bSBarry Smith mat = (Mat_SeqBAIJ *) submats[i]->data; 1191d5b485abSSatish Balay imat_ilen = mat->ilen; 1192d5b485abSSatish Balay imat_j = mat->j; 1193d5b485abSSatish Balay imat_i = mat->i; 1194d5b485abSSatish Balay imat_a = mat->a; 1195d5b485abSSatish Balay cmap_i = cmap[i]; 1196d5b485abSSatish Balay rmap_i = rmap[i]; 1197d5b485abSSatish Balay irow_i = irow[i]; 1198d5b485abSSatish Balay jmax = nrow[i]; 1199d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1200d5b485abSSatish Balay row = irow_i[j]; 1201d5b485abSSatish Balay proc = rtable[row]; 1202d5b485abSSatish Balay if (proc == rank) { 120336f4e84dSSatish Balay row = row - rstart; 120436f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 120536f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 120636f4e84dSSatish Balay cworkA = a_j + a_i[row]; 120736f4e84dSSatish Balay cworkB = b_j + b_i[row]; 120836f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 120936f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 121036f4e84dSSatish Balay 121136f4e84dSSatish Balay row = rmap_i[row + rstart]; 1212df36731bSBarry Smith mat_i = imat_i[row]; 121336f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1214d5b485abSSatish Balay mat_j = imat_j + mat_i; 121536f4e84dSSatish Balay ilen_row = imat_ilen[row]; 121636f4e84dSSatish Balay 121736f4e84dSSatish Balay /* load the column indices for this row into cols*/ 121836f4e84dSSatish Balay for (l=0; l<nzB; l++) { 121936f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 122036f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1221df36731bSBarry Smith *mat_j++ = tcol - 1; 122236f4e84dSSatish Balay PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 1223d5b485abSSatish Balay ilen_row++; 1224d5b485abSSatish Balay } 1225d5b485abSSatish Balay } 122636f4e84dSSatish Balay else break; 122736f4e84dSSatish Balay } 122836f4e84dSSatish Balay imark = l; 122936f4e84dSSatish Balay for (l=0; l<nzA; l++) { 123036f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 123136f4e84dSSatish Balay *mat_j++ = tcol - 1; 123236f4e84dSSatish Balay PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 123336f4e84dSSatish Balay ilen_row++; 123436f4e84dSSatish Balay } 123536f4e84dSSatish Balay } 123636f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 123736f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 123836f4e84dSSatish Balay *mat_j++ = tcol - 1; 123936f4e84dSSatish Balay PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 124036f4e84dSSatish Balay ilen_row++; 124136f4e84dSSatish Balay } 124236f4e84dSSatish Balay } 1243d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1244d5b485abSSatish Balay } 1245d5b485abSSatish Balay } 124636f4e84dSSatish Balay 1247d5b485abSSatish Balay } 1248d5b485abSSatish Balay } 1249d5b485abSSatish Balay 1250d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1251d5b485abSSatish Balay { 1252d5b485abSSatish Balay int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1253d5b485abSSatish Balay int *imat_j,*imat_i; 1254d5b485abSSatish Balay Scalar *imat_a,*rbuf4_i; 1255d5b485abSSatish Balay 1256d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1257d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits4, &i, r_status4+tmp2); 1258d5b485abSSatish Balay index = pa[i]; 1259d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1260d5b485abSSatish Balay jmax = sbuf1_i[0]; 1261d5b485abSSatish Balay ct1 = 2*jmax + 1; 1262d5b485abSSatish Balay ct2 = 0; 1263d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1264d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1265d5b485abSSatish Balay rbuf4_i = rbuf4[i]; 1266d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1267d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1268d5b485abSSatish Balay rmap_i = rmap[is_no]; 1269d5b485abSSatish Balay cmap_i = cmap[is_no]; 1270df36731bSBarry Smith mat = (Mat_SeqBAIJ *) submats[is_no]->data; 1271d5b485abSSatish Balay imat_ilen = mat->ilen; 1272d5b485abSSatish Balay imat_j = mat->j; 1273d5b485abSSatish Balay imat_i = mat->i; 1274d5b485abSSatish Balay imat_a = mat->a; 1275d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1276d5b485abSSatish Balay for (k=0; k<max1; k++, ct1++) { 1277d5b485abSSatish Balay row = sbuf1_i[ct1]; 1278d5b485abSSatish Balay row = rmap_i[row]; 1279d5b485abSSatish Balay ilen = imat_ilen[row]; 1280df36731bSBarry Smith mat_i = imat_i[row]; 128136f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1282d5b485abSSatish Balay mat_j = imat_j + mat_i; 1283d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1284d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1285d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1286df36731bSBarry Smith *mat_j++ = tcol - 1; 128736f4e84dSSatish Balay /* *mat_a++ = rbuf4_i[ct2]; */ 128836f4e84dSSatish Balay PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 1289d5b485abSSatish Balay ilen++; 1290d5b485abSSatish Balay } 1291d5b485abSSatish Balay } 1292d5b485abSSatish Balay imat_ilen[row] = ilen; 1293d5b485abSSatish Balay } 1294d5b485abSSatish Balay } 1295d5b485abSSatish Balay } 1296d5b485abSSatish Balay } 1297d5b485abSSatish Balay PetscFree(r_status4); PetscFree(r_waits4); 1298d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits4,s_status4); 1299d5b485abSSatish Balay PetscFree(s_waits4); PetscFree(s_status4); 1300d5b485abSSatish Balay 1301d5b485abSSatish Balay /* Restore the indices */ 1302d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1303d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i], irow+i); CHKERRQ(ierr); 1304d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i], icol+i); CHKERRQ(ierr); 1305d5b485abSSatish Balay } 1306d5b485abSSatish Balay 1307d5b485abSSatish Balay /* Destroy allocated memory */ 1308d5b485abSSatish Balay PetscFree(irow); 1309d5b485abSSatish Balay PetscFree(w1); 1310d5b485abSSatish Balay PetscFree(pa); 1311d5b485abSSatish Balay 1312d5b485abSSatish Balay PetscFree(sbuf1); 1313d5b485abSSatish Balay PetscFree(rbuf2); 1314d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1315d5b485abSSatish Balay PetscFree(sbuf2[i]); 1316d5b485abSSatish Balay } 1317d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1318d5b485abSSatish Balay PetscFree(rbuf3[i]); 1319d5b485abSSatish Balay PetscFree(rbuf4[i]); 1320d5b485abSSatish Balay } 1321d5b485abSSatish Balay 1322d5b485abSSatish Balay PetscFree(sbuf2); 1323d5b485abSSatish Balay PetscFree(rbuf3); 1324d5b485abSSatish Balay PetscFree(rbuf4 ); 1325d5b485abSSatish Balay PetscFree(sbuf_aj[0]); 1326d5b485abSSatish Balay PetscFree(sbuf_aj); 1327d5b485abSSatish Balay PetscFree(sbuf_aa[0]); 1328d5b485abSSatish Balay PetscFree(sbuf_aa); 1329d5b485abSSatish Balay 1330d5b485abSSatish Balay PetscFree(cmap); 1331d5b485abSSatish Balay PetscFree(rmap); 1332d5b485abSSatish Balay PetscFree(lens); 1333d5b485abSSatish Balay 1334d5b485abSSatish Balay for (i=0; i<ismax; i++) { 133536f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 133636f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1337d5b485abSSatish Balay } 133834e6ae18SSatish Balay 133934e6ae18SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag3); CHKERRQ(ierr); 134034e6ae18SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag2); CHKERRQ(ierr); 134134e6ae18SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1); CHKERRQ(ierr); 134234e6ae18SSatish Balay 1343d5b485abSSatish Balay return 0; 1344d5b485abSSatish Balay } 1345