1d5b485abSSatish Balay #ifndef lint 2*5b83ace0SSatish Balay static char vcid[] = "$Id: baijov.c,v 1.5 1996/07/09 19:02:16 balay Exp balay $"; 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 */ 8df36731bSBarry Smith #include "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 1736f4e84dSSatish Balay static int MatCompressIndices_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 18df36731bSBarry Smith { 19df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 20df36731bSBarry Smith int ierr,isz,bs = baij->bs,Nbs,n,i,j,*idx,*nidx,ival; 21df36731bSBarry Smith char *table; 22df36731bSBarry Smith 23df36731bSBarry Smith Nbs = baij->Nbs; 24df36731bSBarry Smith table = (char *) PetscMalloc((Nbs/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 25df36731bSBarry Smith nidx = (int *) PetscMalloc((Nbs+1)*sizeof(int)); CHKPTRQ(nidx); 26df36731bSBarry Smith 27df36731bSBarry Smith for ( i=0; i<imax; i++ ) { 28df36731bSBarry Smith isz = 0; 29df36731bSBarry Smith PetscMemzero(table,(Nbs/BITSPERBYTE +1)*sizeof(char)); 3036f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 3136f4e84dSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 32df36731bSBarry Smith for ( j=0; j<n ; ++j){ 33df36731bSBarry Smith ival = idx[j]/bs; /* convert the indices into block indices */ 34df36731bSBarry Smith if (ival>Nbs) SETERRQ(1,"MatIncreaseOverlap_SeqBAIJ: index greater than mat-dim"); 35df36731bSBarry Smith if(!BT_LOOKUP(table, ival)) { nidx[isz++] = ival;} 36df36731bSBarry Smith } 3736f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 3836f4e84dSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is_out+i)); CHKERRQ(ierr); 39df36731bSBarry Smith } 40df36731bSBarry Smith PetscFree(table); 41df36731bSBarry Smith PetscFree(nidx); 42df36731bSBarry Smith return 0; 43df36731bSBarry Smith } 44df36731bSBarry Smith 4536f4e84dSSatish Balay static int MatExpandIndices_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 46df36731bSBarry Smith { 47df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 48df36731bSBarry Smith int ierr,bs = baij->bs,Nbs,n,i,j,k,*idx,*nidx; 49df36731bSBarry Smith 50df36731bSBarry Smith Nbs = baij->Nbs; 51df36731bSBarry Smith 52df36731bSBarry Smith nidx = (int *) PetscMalloc((Nbs*bs+1)*sizeof(int)); CHKPTRQ(nidx); 53df36731bSBarry Smith 54df36731bSBarry Smith for ( i=0; i<imax; i++ ) { 5536f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 5636f4e84dSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 57df36731bSBarry Smith for ( j=0; j<n ; ++j){ 58df36731bSBarry Smith for ( k=0; k<bs; k++ ) 59df36731bSBarry Smith nidx[j*bs+k] = idx[j]*bs+k; 60df36731bSBarry Smith } 6136f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 6236f4e84dSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, n*bs, nidx, (is_out+i)); CHKERRQ(ierr); 63df36731bSBarry Smith } 64df36731bSBarry Smith PetscFree(nidx); 65df36731bSBarry Smith return 0; 66df36731bSBarry Smith } 67df36731bSBarry Smith 68df36731bSBarry Smith 69df36731bSBarry Smith int MatIncreaseOverlap_MPIBAIJ(Mat C, int imax, IS *is, int ov) 70d5b485abSSatish Balay { 71d5b485abSSatish Balay int i, ierr; 7236f4e84dSSatish Balay IS *is_new; 7336f4e84dSSatish Balay 7436f4e84dSSatish Balay is_new = (IS *)PetscMalloc(imax*sizeof(IS)); CHKPTRQ(is_new); 75df36731bSBarry Smith /* Convert the indices into block format */ 7636f4e84dSSatish Balay ierr = MatCompressIndices_MPIBAIJ(C, imax, is,is_new); CHKERRQ(ierr); 77df36731bSBarry Smith if (ov < 0){ SETERRQ(1," MatIncreaseOverlap_MPIBAIJ:Negative overlap specified\n");} 78d5b485abSSatish Balay for ( i=0; i<ov; ++i ) { 7936f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new); CHKERRQ(ierr); 80d5b485abSSatish Balay } 8136f4e84dSSatish Balay for ( i=0; i<imax; i++ ) ISDestroy(is[i]); 8236f4e84dSSatish Balay ierr = MatExpandIndices_MPIBAIJ(C, imax, is_new,is); CHKERRQ(ierr); 8336f4e84dSSatish Balay for ( i=0; i<imax; i++ ) ISDestroy(is_new[i]); 8436f4e84dSSatish Balay PetscFree(is_new); 85d5b485abSSatish Balay return 0; 86d5b485abSSatish Balay } 87d5b485abSSatish Balay 88d5b485abSSatish Balay /* 89d5b485abSSatish Balay Sample message format: 90d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 91d5b485abSSatish Balay to index sets 1s[1], is[5] 92d5b485abSSatish Balay mesg [0] = 2 ( no of index sets in the mesg) 93d5b485abSSatish Balay ----------- 94d5b485abSSatish Balay mesg [1] = 1 => is[1] 95d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 96d5b485abSSatish Balay ----------- 97d5b485abSSatish Balay mesg [5] = 5 => is[5] 98d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 99d5b485abSSatish Balay ----------- 100d5b485abSSatish Balay mesg [7] 101d5b485abSSatish Balay mesg [n] datas[1] 102d5b485abSSatish Balay ----------- 103d5b485abSSatish Balay mesg[n+1] 104d5b485abSSatish Balay mesg[m] data(is[5]) 105d5b485abSSatish Balay ----------- 106d5b485abSSatish Balay 107d5b485abSSatish Balay Notes: 108d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 109d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 110d5b485abSSatish Balay */ 111df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat C, int imax, IS *is) 112d5b485abSSatish Balay { 113df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 114d5b485abSSatish Balay int **idx, *n, *w1, *w2, *w3, *w4, *rtable,**data,len,*idx_i; 115df36731bSBarry Smith int size,rank,Mbs,i,j,k,ierr,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 116d5b485abSSatish Balay int *ctr,*pa,tag,*tmp,bsz,nrqr,*isz,*isz1,**xdata,bsz1,**rbuf2; 117d5b485abSSatish Balay char **table; 118d5b485abSSatish Balay MPI_Comm comm; 119d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 120d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 121d5b485abSSatish Balay 122d5b485abSSatish Balay comm = C->comm; 123d5b485abSSatish Balay tag = C->tag; 124d5b485abSSatish Balay size = c->size; 125d5b485abSSatish Balay rank = c->rank; 126df36731bSBarry Smith Mbs = c->Mbs; 127d5b485abSSatish Balay 128df36731bSBarry Smith len = (imax+1)*sizeof(int *) + (imax + Mbs)*sizeof(int); 129d5b485abSSatish Balay idx = (int **) PetscMalloc(len); CHKPTRQ(idx); 130d5b485abSSatish Balay n = (int *) (idx + imax); 131d5b485abSSatish Balay rtable = n + imax; 132d5b485abSSatish Balay 133d5b485abSSatish Balay for ( i=0; i<imax; i++ ) { 134d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]); CHKERRQ(ierr); 135d5b485abSSatish Balay ierr = ISGetSize(is[i],&n[i]); CHKERRQ(ierr); 136d5b485abSSatish Balay } 137d5b485abSSatish Balay 138d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 139d5b485abSSatish Balay for ( i=0,j=0; i<size; i++ ) { 140d5b485abSSatish Balay len = c->rowners[i+1]; 141d5b485abSSatish Balay for ( ; j<len; j++ ) { 142d5b485abSSatish Balay rtable[j] = i; 143d5b485abSSatish Balay } 144d5b485abSSatish Balay } 145d5b485abSSatish Balay 146d5b485abSSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 147d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 148d5b485abSSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1);/* mesg size */ 149d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 150d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 151d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 152d5b485abSSatish Balay PetscMemzero(w1,size*3*sizeof(int)); /* initialise work vector*/ 153d5b485abSSatish Balay for ( i=0; i<imax; i++ ) { 154d5b485abSSatish Balay PetscMemzero(w4,size*sizeof(int)); /* initialise work vector*/ 155d5b485abSSatish Balay idx_i = idx[i]; 156d5b485abSSatish Balay len = n[i]; 157d5b485abSSatish Balay for ( j=0; j<len; j++ ) { 158d5b485abSSatish Balay row = idx_i[j]; 159d5b485abSSatish Balay proc = rtable[row]; 160d5b485abSSatish Balay w4[proc]++; 161d5b485abSSatish Balay } 162d5b485abSSatish Balay for ( j=0; j<size; j++ ){ 163d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 164d5b485abSSatish Balay } 165d5b485abSSatish Balay } 166d5b485abSSatish Balay 167d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 168d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 169d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 170d5b485abSSatish Balay w3[rank] = 0; 171d5b485abSSatish Balay for ( i=0; i<size; i++ ) { 172d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 173d5b485abSSatish Balay } 174d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 175d5b485abSSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); 176d5b485abSSatish Balay for ( i=0,j=0; i<size; i++ ) { 177d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 178d5b485abSSatish Balay } 179d5b485abSSatish Balay 180d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 181d5b485abSSatish Balay for ( i=0; i<nrqs; i++ ) { 182d5b485abSSatish Balay j = pa[i]; 183d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 184d5b485abSSatish Balay msz += w1[j]; 185d5b485abSSatish Balay } 186d5b485abSSatish Balay 187d5b485abSSatish Balay 188d5b485abSSatish Balay /* Do a global reduction to determine how many messages to expect*/ 189d5b485abSSatish Balay { 190d5b485abSSatish Balay int *rw1, *rw2; 191d5b485abSSatish Balay rw1 = (int *) PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 192d5b485abSSatish Balay rw2 = rw1+size; 193d5b485abSSatish Balay MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm); 194d5b485abSSatish Balay bsz = rw1[rank]; 195d5b485abSSatish Balay MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm); 196d5b485abSSatish Balay nrqr = rw2[rank]; 197d5b485abSSatish Balay PetscFree(rw1); 198d5b485abSSatish Balay } 199d5b485abSSatish Balay 200d5b485abSSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 201d5b485abSSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 202d5b485abSSatish Balay rbuf = (int**) PetscMalloc(len); CHKPTRQ(rbuf); 203d5b485abSSatish Balay rbuf[0] = (int *) (rbuf + nrqr); 204d5b485abSSatish Balay for ( i=1; i<nrqr; ++i ) rbuf[i] = rbuf[i-1] + bsz; 205d5b485abSSatish Balay 206d5b485abSSatish Balay /* Post the receives */ 207d5b485abSSatish Balay r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 208d5b485abSSatish Balay CHKPTRQ(r_waits1); 209d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ){ 210d5b485abSSatish Balay MPI_Irecv(rbuf[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i); 211d5b485abSSatish Balay } 212d5b485abSSatish Balay 213d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 214d5b485abSSatish Balay len = 2*size*sizeof(int*) + (size+msz)*sizeof(int); 215d5b485abSSatish Balay outdat = (int **)PetscMalloc(len); CHKPTRQ(outdat); 216d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 217d5b485abSSatish Balay PetscMemzero(outdat,2*size*sizeof(int*)); 218d5b485abSSatish Balay tmp = (int *) (outdat + 2*size); 219d5b485abSSatish Balay ctr = tmp + msz; 220d5b485abSSatish Balay 221d5b485abSSatish Balay { 222d5b485abSSatish Balay int *iptr = tmp,ict = 0; 223d5b485abSSatish Balay for ( i=0; i<nrqs; i++ ) { 224d5b485abSSatish Balay j = pa[i]; 225d5b485abSSatish Balay iptr += ict; 226d5b485abSSatish Balay outdat[j] = iptr; 227d5b485abSSatish Balay ict = w1[j]; 228d5b485abSSatish Balay } 229d5b485abSSatish Balay } 230d5b485abSSatish Balay 231d5b485abSSatish Balay /* Form the outgoing messages */ 232d5b485abSSatish Balay /*plug in the headers*/ 233d5b485abSSatish Balay for ( i=0; i<nrqs; i++ ) { 234d5b485abSSatish Balay j = pa[i]; 235d5b485abSSatish Balay outdat[j][0] = 0; 236d5b485abSSatish Balay PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int)); 237d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 238d5b485abSSatish Balay } 239d5b485abSSatish Balay 240d5b485abSSatish Balay /* Memory for doing local proc's work*/ 241d5b485abSSatish Balay { 242d5b485abSSatish Balay int *d_p; 243d5b485abSSatish Balay char *t_p; 244d5b485abSSatish Balay 245d5b485abSSatish Balay len = (imax+1)*(sizeof(char *) + sizeof(int *) + sizeof(int)) + 246df36731bSBarry Smith (Mbs+1)*imax*sizeof(int) + (Mbs/BITSPERBYTE+1)*imax*sizeof(char); 247d5b485abSSatish Balay table = (char **)PetscMalloc(len); CHKPTRQ(table); 248d5b485abSSatish Balay data = (int **)(table + imax); 249d5b485abSSatish Balay data[0] = (int *)(data + imax); 250df36731bSBarry Smith isz = data[0] + (Mbs+1)*imax; 251d5b485abSSatish Balay table[0] = (char *)(isz + imax); 252d5b485abSSatish Balay d_p = data[0]; t_p = table[0]; 253d5b485abSSatish Balay for ( i=1; i<imax; i++ ) { 254df36731bSBarry Smith table[i] = t_p + (Mbs/BITSPERBYTE+1)*i; 255df36731bSBarry Smith data[i] = d_p + (Mbs+1)*i; 256d5b485abSSatish Balay } 257d5b485abSSatish Balay } 258df36731bSBarry Smith PetscMemzero(*table,(Mbs/BITSPERBYTE+1)*imax); 259d5b485abSSatish Balay PetscMemzero(isz,imax*sizeof(int)); 260d5b485abSSatish Balay 261d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 262d5b485abSSatish Balay { 263d5b485abSSatish Balay int n_i,*data_i,isz_i,*outdat_j,ctr_j; 264d5b485abSSatish Balay char *table_i; 265d5b485abSSatish Balay 266d5b485abSSatish Balay for ( i=0; i<imax; i++ ) { 267d5b485abSSatish Balay PetscMemzero(ctr,size*sizeof(int)); 268d5b485abSSatish Balay n_i = n[i]; 269d5b485abSSatish Balay table_i = table[i]; 270d5b485abSSatish Balay idx_i = idx[i]; 271d5b485abSSatish Balay data_i = data[i]; 272d5b485abSSatish Balay isz_i = isz[i]; 273d5b485abSSatish Balay for ( j=0; j<n_i; j++ ) { /* parse the indices of each IS */ 274d5b485abSSatish Balay row = idx_i[j]; 275d5b485abSSatish Balay proc = rtable[row]; 276d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 277d5b485abSSatish Balay ctr[proc]++; 278d5b485abSSatish Balay *ptr[proc] = row; 279d5b485abSSatish Balay ptr[proc]++; 280d5b485abSSatish Balay } 281d5b485abSSatish Balay else { /* Update the local table */ 282d5b485abSSatish Balay if (!BT_LOOKUP(table_i,row)) { data_i[isz_i++] = row;} 283d5b485abSSatish Balay } 284d5b485abSSatish Balay } 285d5b485abSSatish Balay /* Update the headers for the current IS */ 286d5b485abSSatish Balay for ( j=0; j<size; j++ ) { /* Can Optimise this loop by using pa[] */ 287d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 288d5b485abSSatish Balay outdat_j = outdat[j]; 289d5b485abSSatish Balay k = ++outdat_j[0]; 290d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 291d5b485abSSatish Balay outdat_j[2*k-1] = i; 292d5b485abSSatish Balay } 293d5b485abSSatish Balay } 294d5b485abSSatish Balay isz[i] = isz_i; 295d5b485abSSatish Balay } 296d5b485abSSatish Balay } 297d5b485abSSatish Balay 298d5b485abSSatish Balay 299d5b485abSSatish Balay 300d5b485abSSatish Balay /* Now post the sends */ 301d5b485abSSatish Balay s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 302d5b485abSSatish Balay CHKPTRQ(s_waits1); 303d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 304d5b485abSSatish Balay j = pa[i]; 305d5b485abSSatish Balay MPI_Isend(outdat[j], w1[j], MPI_INT, j, tag, comm, s_waits1+i); 306d5b485abSSatish Balay } 307d5b485abSSatish Balay 308d5b485abSSatish Balay /* No longer need the original indices*/ 309d5b485abSSatish Balay for ( i=0; i<imax; ++i ) { 310d5b485abSSatish Balay ierr = ISRestoreIndices(is[i], idx+i); CHKERRQ(ierr); 311d5b485abSSatish Balay } 312d5b485abSSatish Balay PetscFree(idx); 313d5b485abSSatish Balay 314d5b485abSSatish Balay for ( i=0; i<imax; ++i ) { 315d5b485abSSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 316d5b485abSSatish Balay } 317d5b485abSSatish Balay 318d5b485abSSatish Balay /* Do Local work*/ 319df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 320d5b485abSSatish Balay 321d5b485abSSatish Balay /* Receive messages*/ 322d5b485abSSatish Balay { 323d5b485abSSatish Balay int index; 324d5b485abSSatish Balay 325d5b485abSSatish Balay recv_status = (MPI_Status *) PetscMalloc( (nrqr+1)*sizeof(MPI_Status) ); 326d5b485abSSatish Balay CHKPTRQ(recv_status); 327d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ) { 328d5b485abSSatish Balay MPI_Waitany(nrqr, r_waits1, &index, recv_status+i); 329d5b485abSSatish Balay } 330d5b485abSSatish Balay 331d5b485abSSatish Balay s_status = (MPI_Status *) PetscMalloc( (nrqs+1)*sizeof(MPI_Status) ); 332d5b485abSSatish Balay CHKPTRQ(s_status); 333d5b485abSSatish Balay MPI_Waitall(nrqs,s_waits1,s_status); 334d5b485abSSatish Balay } 335d5b485abSSatish Balay 336d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 337d5b485abSSatish Balay PetscFree(outdat); 338d5b485abSSatish Balay PetscFree(w1); 339d5b485abSSatish Balay 340d5b485abSSatish Balay xdata = (int **)PetscMalloc((nrqr+1)*sizeof(int *)); CHKPTRQ(xdata); 341d5b485abSSatish Balay isz1 = (int *)PetscMalloc((nrqr+1)*sizeof(int)); CHKPTRQ(isz1); 342df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 343d5b485abSSatish Balay PetscFree(rbuf); 344d5b485abSSatish Balay 345d5b485abSSatish Balay /* Send the data back*/ 346d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 347d5b485abSSatish Balay { 348d5b485abSSatish Balay int *rw1, *rw2; 349d5b485abSSatish Balay 350d5b485abSSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 351d5b485abSSatish Balay PetscMemzero(rw1,2*size*sizeof(int)); 352d5b485abSSatish Balay rw2 = rw1+size; 353d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ) { 354d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 355d5b485abSSatish Balay rw1[proc] = isz1[i]; 356d5b485abSSatish Balay } 357d5b485abSSatish Balay 358d5b485abSSatish Balay MPI_Allreduce(rw1, rw2, size, MPI_INT, MPI_MAX, comm); 359d5b485abSSatish Balay bsz1 = rw2[rank]; 360d5b485abSSatish Balay PetscFree(rw1); 361d5b485abSSatish Balay } 362d5b485abSSatish Balay 363d5b485abSSatish Balay /* Allocate buffers*/ 364d5b485abSSatish Balay 365d5b485abSSatish Balay /* Allocate memory for recv buffers. Prob none if nrqr = 0 ???? */ 366d5b485abSSatish Balay len = (nrqs+1)*sizeof(int*) + nrqs*bsz1*sizeof(int); 367d5b485abSSatish Balay rbuf2 = (int**) PetscMalloc(len); CHKPTRQ(rbuf2); 368d5b485abSSatish Balay rbuf2[0] = (int *) (rbuf2 + nrqs); 369d5b485abSSatish Balay for ( i=1; i<nrqs; ++i ) rbuf2[i] = rbuf2[i-1] + bsz1; 370d5b485abSSatish Balay 371d5b485abSSatish Balay /* Post the receives */ 372d5b485abSSatish Balay r_waits2 = (MPI_Request *)PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 373d5b485abSSatish Balay CHKPTRQ(r_waits2); 374d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 375d5b485abSSatish Balay MPI_Irecv(rbuf2[i], bsz1, MPI_INT, MPI_ANY_SOURCE, tag, comm, r_waits2+i); 376d5b485abSSatish Balay } 377d5b485abSSatish Balay 378d5b485abSSatish Balay /* Now post the sends */ 379d5b485abSSatish Balay s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 380d5b485abSSatish Balay CHKPTRQ(s_waits2); 381d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ) { 382d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 383d5b485abSSatish Balay MPI_Isend( xdata[i], isz1[i], MPI_INT, j, tag, comm, s_waits2+i); 384d5b485abSSatish Balay } 385d5b485abSSatish Balay 386d5b485abSSatish Balay /* receive work done on other processors*/ 387d5b485abSSatish Balay { 388d5b485abSSatish Balay int index, is_no, ct1, max,*rbuf2_i,isz_i,*data_i,jmax; 389d5b485abSSatish Balay char *table_i; 390d5b485abSSatish Balay MPI_Status *status2; 391d5b485abSSatish Balay 392d5b485abSSatish Balay status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(status2); 393d5b485abSSatish Balay 394d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 395d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits2, &index, status2+i); 396d5b485abSSatish Balay /* Process the message*/ 397d5b485abSSatish Balay rbuf2_i = rbuf2[index]; 398d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 399d5b485abSSatish Balay jmax = rbuf2[index][0]; 400d5b485abSSatish Balay for ( j=1; j<=jmax; j++ ) { 401d5b485abSSatish Balay max = rbuf2_i[2*j]; 402d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 403d5b485abSSatish Balay isz_i = isz[is_no]; 404d5b485abSSatish Balay data_i = data[is_no]; 405d5b485abSSatish Balay table_i = table[is_no]; 406d5b485abSSatish Balay for ( k=0; k<max; k++,ct1++ ) { 407d5b485abSSatish Balay row = rbuf2_i[ct1]; 408d5b485abSSatish Balay if (!BT_LOOKUP(table_i,row)) { data_i[isz_i++] = row;} 409d5b485abSSatish Balay } 410d5b485abSSatish Balay isz[is_no] = isz_i; 411d5b485abSSatish Balay } 412d5b485abSSatish Balay } 413d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits2,status2); 414d5b485abSSatish Balay PetscFree(status2); 415d5b485abSSatish Balay } 416d5b485abSSatish Balay 417d5b485abSSatish Balay for ( i=0; i<imax; ++i ) { 418d5b485abSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz[i], data[i], is+i); CHKERRQ(ierr); 419d5b485abSSatish Balay } 420d5b485abSSatish Balay 421d5b485abSSatish Balay PetscFree(pa); 422d5b485abSSatish Balay PetscFree(rbuf2); 423d5b485abSSatish Balay PetscFree(s_waits1); 424d5b485abSSatish Balay PetscFree(r_waits1); 425d5b485abSSatish Balay PetscFree(s_waits2); 426d5b485abSSatish Balay PetscFree(r_waits2); 427d5b485abSSatish Balay PetscFree(table); 428d5b485abSSatish Balay PetscFree(s_status); 429d5b485abSSatish Balay PetscFree(recv_status); 430d5b485abSSatish Balay PetscFree(xdata[0]); 431d5b485abSSatish Balay PetscFree(xdata); 432d5b485abSSatish Balay PetscFree(isz1); 433d5b485abSSatish Balay return 0; 434d5b485abSSatish Balay } 435d5b485abSSatish Balay 436d5b485abSSatish Balay /* 437df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 438d5b485abSSatish Balay the work on the local processor. 439d5b485abSSatish Balay 440d5b485abSSatish Balay Inputs: 441df36731bSBarry Smith C - MAT_MPIBAIJ; 442d5b485abSSatish Balay imax - total no of index sets processed at a time; 443df36731bSBarry Smith table - an array of char - size = Mbs bits. 444d5b485abSSatish Balay 445d5b485abSSatish Balay Output: 446d5b485abSSatish Balay isz - array containing the count of the solution elements correspondign 447d5b485abSSatish Balay to each index set; 448d5b485abSSatish Balay data - pointer to the solutions 449d5b485abSSatish Balay */ 450df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,char **table,int *isz, 451d5b485abSSatish Balay int **data) 452d5b485abSSatish Balay { 453df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 454d5b485abSSatish Balay Mat A = c->A, B = c->B; 455df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 456df36731bSBarry Smith int start, end, val, max, rstart,cstart,*ai, *aj; 457d5b485abSSatish Balay int *bi, *bj, *garray, i, j, k, row,*data_i,isz_i; 458d5b485abSSatish Balay char *table_i; 459d5b485abSSatish Balay 460d5b485abSSatish Balay rstart = c->rstart; 461d5b485abSSatish Balay cstart = c->cstart; 462d5b485abSSatish Balay ai = a->i; 463df36731bSBarry Smith aj = a->j; 464d5b485abSSatish Balay bi = b->i; 465df36731bSBarry Smith bj = b->j; 466d5b485abSSatish Balay garray = c->garray; 467d5b485abSSatish Balay 468d5b485abSSatish Balay 469d5b485abSSatish Balay for ( i=0; i<imax; i++ ) { 470d5b485abSSatish Balay data_i = data[i]; 471d5b485abSSatish Balay table_i = table[i]; 472d5b485abSSatish Balay isz_i = isz[i]; 473d5b485abSSatish Balay for ( j=0, max=isz[i]; j<max; j++ ) { 474d5b485abSSatish Balay row = data_i[j] - rstart; 475d5b485abSSatish Balay start = ai[row]; 476d5b485abSSatish Balay end = ai[row+1]; 477d5b485abSSatish Balay for ( k=start; k<end; k++ ) { /* Amat */ 478df36731bSBarry Smith val = aj[k] + cstart; 479d5b485abSSatish Balay if (!BT_LOOKUP(table_i,val)) { data_i[isz_i++] = val;} 480d5b485abSSatish Balay } 481d5b485abSSatish Balay start = bi[row]; 482d5b485abSSatish Balay end = bi[row+1]; 483d5b485abSSatish Balay for ( k=start; k<end; k++ ) { /* Bmat */ 484df36731bSBarry Smith val = garray[bj[k]]; 485d5b485abSSatish Balay if (!BT_LOOKUP(table_i,val)) { data_i[isz_i++] = val;} 486d5b485abSSatish Balay } 487d5b485abSSatish Balay } 488d5b485abSSatish Balay isz[i] = isz_i; 489d5b485abSSatish Balay } 490d5b485abSSatish Balay return 0; 491d5b485abSSatish Balay } 492d5b485abSSatish Balay /* 493df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 494d5b485abSSatish Balay and return the output 495d5b485abSSatish Balay 496d5b485abSSatish Balay Input: 497d5b485abSSatish Balay C - the matrix 498d5b485abSSatish Balay nrqr - no of messages being processed. 499d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 500d5b485abSSatish Balay 501d5b485abSSatish Balay Output: 502d5b485abSSatish Balay xdata - array of messages to be sent back 503d5b485abSSatish Balay isz1 - size of each message 504d5b485abSSatish Balay 505d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 506d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 507d5b485abSSatish Balay rather then all previous rows as it is now where a single large chunck of 508d5b485abSSatish Balay memory is used. 509d5b485abSSatish Balay 510d5b485abSSatish Balay */ 511df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf, 512d5b485abSSatish Balay int **xdata, int * isz1) 513d5b485abSSatish Balay { 514df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 515d5b485abSSatish Balay Mat A = c->A, B = c->B; 516df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 517df36731bSBarry Smith int rstart,cstart,*ai, *aj, *bi, *bj, *garray, i, j, k; 518d5b485abSSatish Balay int row,total_sz,ct, ct1, ct2, ct3,mem_estimate, oct2, l, start, end; 519df36731bSBarry Smith int val, max1, max2, rank, Mbs, no_malloc =0, *tmp, new_estimate, ctr; 520d5b485abSSatish Balay int *rbuf_i,kmax,rbuf_0; 521d5b485abSSatish Balay char *xtable; 522d5b485abSSatish Balay 523d5b485abSSatish Balay rank = c->rank; 524df36731bSBarry Smith Mbs = c->Mbs; 525d5b485abSSatish Balay rstart = c->rstart; 526d5b485abSSatish Balay cstart = c->cstart; 527d5b485abSSatish Balay ai = a->i; 528df36731bSBarry Smith aj = a->j; 529d5b485abSSatish Balay bi = b->i; 530df36731bSBarry Smith bj = b->j; 531d5b485abSSatish Balay garray = c->garray; 532d5b485abSSatish Balay 533d5b485abSSatish Balay 534d5b485abSSatish Balay for ( i=0,ct=0,total_sz=0; i<nrqr; ++i ) { 535d5b485abSSatish Balay rbuf_i = rbuf[i]; 536d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 537d5b485abSSatish Balay ct += rbuf_0; 538d5b485abSSatish Balay for ( j=1; j<=rbuf_0; j++ ) { total_sz += rbuf_i[2*j]; } 539d5b485abSSatish Balay } 540d5b485abSSatish Balay 541df36731bSBarry Smith max1 = ct*(a->nz +b->nz)/c->Mbs; 542d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 543d5b485abSSatish Balay xdata[0] = (int *)PetscMalloc(mem_estimate*sizeof(int)); CHKPTRQ(xdata[0]); 544d5b485abSSatish Balay ++no_malloc; 545df36731bSBarry Smith xtable = (char *)PetscMalloc((Mbs/BITSPERBYTE+1)*sizeof(char)); CHKPTRQ(xtable); 546d5b485abSSatish Balay PetscMemzero(isz1,nrqr*sizeof(int)); 547d5b485abSSatish Balay 548d5b485abSSatish Balay ct3 = 0; 549d5b485abSSatish Balay for ( i=0; i<nrqr; i++ ) { /* for easch mesg from proc i */ 550d5b485abSSatish Balay rbuf_i = rbuf[i]; 551d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 552d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 553d5b485abSSatish Balay ct2 = ct1; 554d5b485abSSatish Balay ct3 += ct1; 555d5b485abSSatish Balay for ( j=1; j<=rbuf_0; j++ ) { /* for each IS from proc i*/ 556df36731bSBarry Smith PetscMemzero(xtable,(Mbs/BITSPERBYTE+1)*sizeof(char)); 557d5b485abSSatish Balay oct2 = ct2; 558d5b485abSSatish Balay kmax = rbuf_i[2*j]; 559d5b485abSSatish Balay for ( k=0; k<kmax; k++, ct1++ ) { 560d5b485abSSatish Balay row = rbuf_i[ct1]; 561d5b485abSSatish Balay if (!BT_LOOKUP(xtable,row)) { 562d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 563d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 564d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 565d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 566d5b485abSSatish Balay PetscFree(xdata[0]); 567d5b485abSSatish Balay xdata[0] = tmp; 568d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 569d5b485abSSatish Balay for ( ctr=1; ctr<=i; ctr++ ) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 570d5b485abSSatish Balay } 571d5b485abSSatish Balay xdata[i][ct2++] = row; 572d5b485abSSatish Balay ct3++; 573d5b485abSSatish Balay } 574d5b485abSSatish Balay } 575d5b485abSSatish Balay for ( k=oct2,max2=ct2; k<max2; k++ ) { 576d5b485abSSatish Balay row = xdata[i][k] - rstart; 577d5b485abSSatish Balay start = ai[row]; 578d5b485abSSatish Balay end = ai[row+1]; 579d5b485abSSatish Balay for ( l=start; l<end; l++ ) { 580df36731bSBarry Smith val = aj[l] + cstart; 581d5b485abSSatish Balay if (!BT_LOOKUP(xtable,val)) { 582d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 583d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 584d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 585d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 586d5b485abSSatish Balay PetscFree(xdata[0]); 587d5b485abSSatish Balay xdata[0] = tmp; 588d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 589d5b485abSSatish Balay for ( ctr=1; ctr<=i; ctr++ ) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 590d5b485abSSatish Balay } 591d5b485abSSatish Balay xdata[i][ct2++] = val; 592d5b485abSSatish Balay ct3++; 593d5b485abSSatish Balay } 594d5b485abSSatish Balay } 595d5b485abSSatish Balay start = bi[row]; 596d5b485abSSatish Balay end = bi[row+1]; 597d5b485abSSatish Balay for ( l=start; l<end; l++ ) { 598df36731bSBarry Smith val = garray[bj[l]]; 599d5b485abSSatish Balay if (!BT_LOOKUP(xtable,val)) { 600d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 601d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 602d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 603d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 604d5b485abSSatish Balay PetscFree(xdata[0]); 605d5b485abSSatish Balay xdata[0] = tmp; 606d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 607d5b485abSSatish Balay for ( ctr =1; ctr <=i; ctr++ ) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 608d5b485abSSatish Balay } 609d5b485abSSatish Balay xdata[i][ct2++] = val; 610d5b485abSSatish Balay ct3++; 611d5b485abSSatish Balay } 612d5b485abSSatish Balay } 613d5b485abSSatish Balay } 614d5b485abSSatish Balay /* Update the header*/ 615d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 616d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 617d5b485abSSatish Balay } 618d5b485abSSatish Balay xdata[i][0] = rbuf_0; 619d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 620d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 621d5b485abSSatish Balay } 622d5b485abSSatish Balay PetscFree(xtable); 623df36731bSBarry Smith PLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d bytes, no of mallocs = %d\n",rank,mem_estimate, ct3,no_malloc); 624d5b485abSSatish Balay return 0; 625d5b485abSSatish Balay } 626d5b485abSSatish Balay 627a2feddc5SSatish Balay static int MatGetSubMatrices_MPIBAIJ_local(Mat,int,IS *,IS *,MatGetSubMatrixCall,Mat **); 628a2feddc5SSatish Balay 629df36731bSBarry Smith int MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,IS *isrow,IS *iscol, 630d5b485abSSatish Balay MatGetSubMatrixCall scall,Mat **submat) 631d5b485abSSatish Balay { 63236f4e84dSSatish Balay int i,ierr; 63336f4e84dSSatish Balay IS *isrow_new,*iscol_new; 634a2feddc5SSatish Balay 635a2feddc5SSatish Balay /* The compression and expansion should be avoided. Dos'nt point 636a2feddc5SSatish Balay out errors might change the indices hence buggey */ 637a2feddc5SSatish Balay 63836f4e84dSSatish Balay isrow_new = (IS *)PetscMalloc(2*ismax*sizeof(IS)); CHKPTRQ(isrow_new); 63936f4e84dSSatish Balay iscol_new = isrow_new + ismax; 64036f4e84dSSatish Balay ierr = MatCompressIndices_MPIBAIJ(C, ismax, isrow,isrow_new); CHKERRQ(ierr); 64136f4e84dSSatish Balay ierr = MatCompressIndices_MPIBAIJ(C, ismax, iscol,iscol_new); CHKERRQ(ierr); 64236f4e84dSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,ismax,isrow_new,iscol_new,scall,submat); CHKERRQ(ierr); 64336f4e84dSSatish Balay 64436f4e84dSSatish Balay for ( i=0; i<ismax; i++) { 64536f4e84dSSatish Balay ISDestroy(isrow_new[i]); 64636f4e84dSSatish Balay ISDestroy(iscol_new[i]); 64736f4e84dSSatish Balay } 64836f4e84dSSatish Balay PetscFree(isrow_new); 649a2feddc5SSatish Balay return 0; 650a2feddc5SSatish Balay } 651a2feddc5SSatish Balay 652a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 653a2feddc5SSatish Balay static int MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,IS *isrow,IS *iscol, 654a2feddc5SSatish Balay MatGetSubMatrixCall scall,Mat **submat) 655a2feddc5SSatish Balay { 656df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 657d5b485abSSatish Balay Mat A = c->A,*submats = *submat; 658df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data, *b = (Mat_SeqBAIJ*)c->B->data, *mat; 659d5b485abSSatish Balay int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size; 660a2feddc5SSatish Balay int **sbuf1,**sbuf2, rank, Mbs,i,j,k,l,ct1,ct2,ierr, **rbuf1,row,proc; 661d5b485abSSatish Balay int nrqs, msz, **ptr,index,*req_size,*ctr,*pa,tag,*tmp,tcol,bsz,nrqr; 662df36731bSBarry Smith int **rbuf3,*req_source,**sbuf_aj, **rbuf2, max1,max2,**rmap; 663d5b485abSSatish Balay int **cmap,**lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 664d5b485abSSatish Balay int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*cmap_i,*lens_i; 66536f4e84dSSatish Balay int *rmap_i,bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA, *cworkB; 66636f4e84dSSatish Balay int cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 66736f4e84dSSatish Balay int *bmap = c->garray,ctmp,rstart=c->rstart; 668d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 669d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 670d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 671d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 672d5b485abSSatish Balay MPI_Comm comm; 67336f4e84dSSatish Balay Scalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 67436f4e84dSSatish Balay Scalar *a_a=a->a,*b_a=b->a; 675d5b485abSSatish Balay 676d5b485abSSatish Balay comm = C->comm; 677d5b485abSSatish Balay tag = C->tag; 678d5b485abSSatish Balay size = c->size; 679d5b485abSSatish Balay rank = c->rank; 680a2feddc5SSatish Balay Mbs = c->Mbs; 681d5b485abSSatish Balay 682d5b485abSSatish Balay /* Check if the col indices are sorted */ 683d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 684d5b485abSSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j); 685df36731bSBarry Smith if (!j) SETERRQ(1,"MatGetSubmatrices_MPIBAIJ:IS is not sorted"); 686d5b485abSSatish Balay } 687d5b485abSSatish Balay 688a2feddc5SSatish Balay len = (2*ismax+1)*(sizeof(int *) + sizeof(int)) + (Mbs+1)*sizeof(int); 689d5b485abSSatish Balay irow = (int **)PetscMalloc(len); CHKPTRQ(irow); 690d5b485abSSatish Balay icol = irow + ismax; 691d5b485abSSatish Balay nrow = (int *) (icol + ismax); 692d5b485abSSatish Balay ncol = nrow + ismax; 693d5b485abSSatish Balay rtable = ncol + ismax; 694d5b485abSSatish Balay 695d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 696d5b485abSSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]); CHKERRQ(ierr); 697d5b485abSSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]); CHKERRQ(ierr); 698d5b485abSSatish Balay ierr = ISGetSize(isrow[i],&nrow[i]); CHKERRQ(ierr); 699d5b485abSSatish Balay ierr = ISGetSize(iscol[i],&ncol[i]); CHKERRQ(ierr); 700d5b485abSSatish Balay } 701d5b485abSSatish Balay 702d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 703d5b485abSSatish Balay for ( i=0,j=0; i<size; i++ ) { 704d5b485abSSatish Balay jmax = c->rowners[i+1]; 705d5b485abSSatish Balay for ( ; j<jmax; j++ ) { 706d5b485abSSatish Balay rtable[j] = i; 707d5b485abSSatish Balay } 708d5b485abSSatish Balay } 709d5b485abSSatish Balay 710d5b485abSSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 711d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 712d5b485abSSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */ 713d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 714d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 715d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 716d5b485abSSatish Balay PetscMemzero(w1,size*3*sizeof(int)); /* initialise work vector*/ 717d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 718d5b485abSSatish Balay PetscMemzero(w4,size*sizeof(int)); /* initialise work vector*/ 719d5b485abSSatish Balay jmax = nrow[i]; 720d5b485abSSatish Balay irow_i = irow[i]; 721d5b485abSSatish Balay for ( j=0; j<jmax; j++ ) { 722d5b485abSSatish Balay row = irow_i[j]; 723d5b485abSSatish Balay proc = rtable[row]; 724d5b485abSSatish Balay w4[proc]++; 725d5b485abSSatish Balay } 726d5b485abSSatish Balay for ( j=0; j<size; j++ ) { 727d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 728d5b485abSSatish Balay } 729d5b485abSSatish Balay } 730d5b485abSSatish Balay 731d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 732d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 733d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 734d5b485abSSatish Balay w3[rank] = 0; 735d5b485abSSatish Balay for ( i=0; i<size; i++ ) { 736d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 737d5b485abSSatish Balay } 738d5b485abSSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/ 739d5b485abSSatish Balay for ( i=0, j=0; i<size; i++ ) { 740d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 741d5b485abSSatish Balay } 742d5b485abSSatish Balay 743d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 744d5b485abSSatish Balay for ( i=0; i<nrqs; i++ ) { 745d5b485abSSatish Balay j = pa[i]; 746d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 747d5b485abSSatish Balay msz += w1[j]; 748d5b485abSSatish Balay } 749d5b485abSSatish Balay /* Do a global reduction to determine how many messages to expect*/ 750d5b485abSSatish Balay { 751d5b485abSSatish Balay int *rw1, *rw2; 752d5b485abSSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 753d5b485abSSatish Balay rw2 = rw1+size; 754d5b485abSSatish Balay MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm); 755d5b485abSSatish Balay bsz = rw1[rank]; 756d5b485abSSatish Balay MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm); 757d5b485abSSatish Balay nrqr = rw2[rank]; 758d5b485abSSatish Balay PetscFree(rw1); 759d5b485abSSatish Balay } 760d5b485abSSatish Balay 761d5b485abSSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 762d5b485abSSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 763d5b485abSSatish Balay rbuf1 = (int**) PetscMalloc(len); CHKPTRQ(rbuf1); 764d5b485abSSatish Balay rbuf1[0] = (int *) (rbuf1 + nrqr); 765d5b485abSSatish Balay for ( i=1; i<nrqr; ++i ) rbuf1[i] = rbuf1[i-1] + bsz; 766d5b485abSSatish Balay 767d5b485abSSatish Balay /* Post the receives */ 768d5b485abSSatish Balay r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 769d5b485abSSatish Balay CHKPTRQ(r_waits1); 770d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ) { 771d5b485abSSatish Balay MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i); 772d5b485abSSatish Balay } 773d5b485abSSatish Balay 774d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 775d5b485abSSatish Balay len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 776d5b485abSSatish Balay sbuf1 = (int **)PetscMalloc(len); CHKPTRQ(sbuf1); 777d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 778d5b485abSSatish Balay PetscMemzero(sbuf1,2*size*sizeof(int*)); 779d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 780d5b485abSSatish Balay tmp = (int *) (ptr + size); 781d5b485abSSatish Balay ctr = tmp + 2*msz; 782d5b485abSSatish Balay 783d5b485abSSatish Balay { 78436f4e84dSSatish Balay 785d5b485abSSatish Balay int *iptr = tmp,ict = 0; 786d5b485abSSatish Balay for ( i=0; i<nrqs; i++ ) { 787d5b485abSSatish Balay j = pa[i]; 788d5b485abSSatish Balay iptr += ict; 789d5b485abSSatish Balay sbuf1[j] = iptr; 790d5b485abSSatish Balay ict = w1[j]; 791d5b485abSSatish Balay } 792d5b485abSSatish Balay } 793d5b485abSSatish Balay 794d5b485abSSatish Balay /* Form the outgoing messages */ 795d5b485abSSatish Balay /* Initialise the header space */ 796d5b485abSSatish Balay for ( i=0; i<nrqs; i++ ) { 797d5b485abSSatish Balay j = pa[i]; 798d5b485abSSatish Balay sbuf1[j][0] = 0; 799d5b485abSSatish Balay PetscMemzero(sbuf1[j]+1, 2*w3[j]*sizeof(int)); 800d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 801d5b485abSSatish Balay } 802d5b485abSSatish Balay 803d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 804d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 805d5b485abSSatish Balay PetscMemzero(ctr,size*sizeof(int)); 806d5b485abSSatish Balay irow_i = irow[i]; 807d5b485abSSatish Balay jmax = nrow[i]; 808d5b485abSSatish Balay for ( j=0; j<jmax; j++ ) { /* parse the indices of each IS */ 809d5b485abSSatish Balay row = irow_i[j]; 810d5b485abSSatish Balay proc = rtable[row]; 811d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 812d5b485abSSatish Balay ctr[proc]++; 813d5b485abSSatish Balay *ptr[proc] = row; 814d5b485abSSatish Balay ptr[proc]++; 815d5b485abSSatish Balay } 816d5b485abSSatish Balay } 817d5b485abSSatish Balay /* Update the headers for the current IS */ 818d5b485abSSatish Balay for ( j=0; j<size; j++ ) { /* Can Optimise this loop too */ 819d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 820d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 821d5b485abSSatish Balay k = ++sbuf1_j[0]; 822d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 823d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 824d5b485abSSatish Balay } 825d5b485abSSatish Balay } 826d5b485abSSatish Balay } 827d5b485abSSatish Balay 828d5b485abSSatish Balay /* Now post the sends */ 829d5b485abSSatish Balay s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 830d5b485abSSatish Balay CHKPTRQ(s_waits1); 831d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 832d5b485abSSatish Balay j = pa[i]; 833d5b485abSSatish Balay /* printf("[%d] Send Req to %d: size %d \n", rank,j, w1[j]); */ 834d5b485abSSatish Balay MPI_Isend( sbuf1[j], w1[j], MPI_INT, j, tag, comm, s_waits1+i); 835d5b485abSSatish Balay } 836d5b485abSSatish Balay 837d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 838d5b485abSSatish Balay r_waits2 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 839d5b485abSSatish Balay CHKPTRQ(r_waits2); 840d5b485abSSatish Balay rbuf2 = (int**)PetscMalloc((nrqs+1)*sizeof(int *));CHKPTRQ(rbuf2); 841d5b485abSSatish Balay rbuf2[0] = tmp + msz; 842d5b485abSSatish Balay for ( i=1; i<nrqs; ++i ) { 843d5b485abSSatish Balay j = pa[i]; 844d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 845d5b485abSSatish Balay } 846d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 847d5b485abSSatish Balay j = pa[i]; 848d5b485abSSatish Balay MPI_Irecv( rbuf2[i], w1[j], MPI_INT, j, tag+1, comm, r_waits2+i); 849d5b485abSSatish Balay } 850d5b485abSSatish Balay 851d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 852d5b485abSSatish Balay 853d5b485abSSatish Balay 854d5b485abSSatish Balay /* Receive messages*/ 855d5b485abSSatish Balay s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 856d5b485abSSatish Balay CHKPTRQ(s_waits2); 857d5b485abSSatish Balay r_status1 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 858d5b485abSSatish Balay CHKPTRQ(r_status1); 859d5b485abSSatish Balay len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*); 860d5b485abSSatish Balay sbuf2 = (int**) PetscMalloc(len);CHKPTRQ(sbuf2); 861d5b485abSSatish Balay req_size = (int *) (sbuf2 + nrqr); 862d5b485abSSatish Balay req_source = req_size + nrqr; 863d5b485abSSatish Balay 864d5b485abSSatish Balay { 865df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*) c->A->data, *sB = (Mat_SeqBAIJ*) c->B->data; 86636f4e84dSSatish Balay int *sAi = sA->i, *sBi = sB->i, id; 867d5b485abSSatish Balay int *sbuf2_i; 868d5b485abSSatish Balay 869d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ) { 870d5b485abSSatish Balay MPI_Waitany(nrqr, r_waits1, &index, r_status1+i); 871d5b485abSSatish Balay req_size[index] = 0; 872d5b485abSSatish Balay rbuf1_i = rbuf1[index]; 873d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 874d5b485abSSatish Balay MPI_Get_count(r_status1+i,MPI_INT, &end); 875d5b485abSSatish Balay sbuf2[index] = (int *)PetscMalloc(end*sizeof(int));CHKPTRQ(sbuf2[index]); 876d5b485abSSatish Balay sbuf2_i = sbuf2[index]; 877d5b485abSSatish Balay for ( j=start; j<end; j++ ) { 878d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 879d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 880d5b485abSSatish Balay sbuf2_i[j] = ncols; 881d5b485abSSatish Balay req_size[index] += ncols; 882d5b485abSSatish Balay } 883d5b485abSSatish Balay req_source[index] = r_status1[i].MPI_SOURCE; 884d5b485abSSatish Balay /* form the header */ 885d5b485abSSatish Balay sbuf2_i[0] = req_size[index]; 886d5b485abSSatish Balay for ( j=1; j<start; j++ ) { sbuf2_i[j] = rbuf1_i[j]; } 887d5b485abSSatish Balay MPI_Isend(sbuf2_i,end,MPI_INT,req_source[index],tag+1,comm,s_waits2+i); 888d5b485abSSatish Balay } 889d5b485abSSatish Balay } 890d5b485abSSatish Balay PetscFree(r_status1); PetscFree(r_waits1); 891d5b485abSSatish Balay 892d5b485abSSatish Balay /* recv buffer sizes */ 893d5b485abSSatish Balay /* Receive messages*/ 894d5b485abSSatish Balay 895d5b485abSSatish Balay rbuf3 = (int**)PetscMalloc((nrqs+1)*sizeof(int*)); CHKPTRQ(rbuf3); 896d5b485abSSatish Balay rbuf4 = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar*));CHKPTRQ(rbuf4); 897d5b485abSSatish Balay r_waits3 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 898d5b485abSSatish Balay CHKPTRQ(r_waits3); 899d5b485abSSatish Balay r_waits4 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request)); 900d5b485abSSatish Balay CHKPTRQ(r_waits4); 901d5b485abSSatish Balay r_status2 = (MPI_Status *) PetscMalloc( (nrqs+1)*sizeof(MPI_Status) ); 902d5b485abSSatish Balay CHKPTRQ(r_status2); 903d5b485abSSatish Balay 904d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 905d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits2, &index, r_status2+i); 906d5b485abSSatish Balay rbuf3[index] = (int *)PetscMalloc(rbuf2[index][0]*sizeof(int)); 907d5b485abSSatish Balay CHKPTRQ(rbuf3[index]); 908a2feddc5SSatish Balay rbuf4[index] = (Scalar *)PetscMalloc(rbuf2[index][0]*bs2*sizeof(Scalar)); 909d5b485abSSatish Balay CHKPTRQ(rbuf4[index]); 910d5b485abSSatish Balay MPI_Irecv(rbuf3[index],rbuf2[index][0], MPI_INT, 911d5b485abSSatish Balay r_status2[i].MPI_SOURCE, tag+2, comm, r_waits3+index); 912a2feddc5SSatish Balay MPI_Irecv(rbuf4[index],rbuf2[index][0]*bs2, MPIU_SCALAR, 913d5b485abSSatish Balay r_status2[i].MPI_SOURCE, tag+3, comm, r_waits4+index); 914d5b485abSSatish Balay } 915d5b485abSSatish Balay PetscFree(r_status2); PetscFree(r_waits2); 916d5b485abSSatish Balay 917d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 918d5b485abSSatish Balay s_status1 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status)); 919d5b485abSSatish Balay CHKPTRQ(s_status1); 920d5b485abSSatish Balay s_status2 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 921d5b485abSSatish Balay CHKPTRQ(s_status2); 922d5b485abSSatish Balay 923d5b485abSSatish Balay MPI_Waitall(nrqs,s_waits1,s_status1); 924d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits2,s_status2); 925d5b485abSSatish Balay PetscFree(s_status1); PetscFree(s_status2); 926d5b485abSSatish Balay PetscFree(s_waits1); PetscFree(s_waits2); 927d5b485abSSatish Balay 928d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 929d5b485abSSatish Balay sbuf_aj = (int **)PetscMalloc((nrqr+1)*sizeof(int *));CHKPTRQ(sbuf_aj); 930d5b485abSSatish Balay for ( i=0,j=0; i<nrqr; i++ ) j += req_size[i]; 931d5b485abSSatish Balay sbuf_aj[0] = (int*) PetscMalloc((j+1)*sizeof(int)); CHKPTRQ(sbuf_aj[0]); 932d5b485abSSatish Balay for ( i=1; i<nrqr; i++ ) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 933d5b485abSSatish Balay 934d5b485abSSatish Balay s_waits3 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 935d5b485abSSatish Balay CHKPTRQ(s_waits3); 936d5b485abSSatish Balay { 937d5b485abSSatish Balay for ( i=0; i<nrqr; i++ ) { 938d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 939d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 940d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 941d5b485abSSatish Balay ct2 = 0; 942d5b485abSSatish Balay for ( j=1,max1=rbuf1_i[0]; j<=max1; j++ ) { 943d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 944d5b485abSSatish Balay for ( k=0; k<kmax; k++,ct1++ ) { 945d5b485abSSatish Balay row = rbuf1_i[ct1] - cstart; 946d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 947d5b485abSSatish Balay ncols = nzA + nzB; 948d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 949d5b485abSSatish Balay 950d5b485abSSatish Balay /* load the column indices for this row into cols*/ 951d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 952d5b485abSSatish Balay for ( l=0; l<nzB; l++ ) { 953d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 954d5b485abSSatish Balay else break; 955d5b485abSSatish Balay } 956d5b485abSSatish Balay imark = l; 957d5b485abSSatish Balay for ( l=0; l<nzA; l++ ) cols[imark+l] = cstart + cworkA[l]; 958d5b485abSSatish Balay for ( l=imark; l<nzB; l++ ) cols[nzA+l] = bmap[cworkB[l]]; 959d5b485abSSatish Balay ct2 += ncols; 960d5b485abSSatish Balay } 961d5b485abSSatish Balay } 962d5b485abSSatish Balay MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag+2,comm,s_waits3+i); 963d5b485abSSatish Balay } 964d5b485abSSatish Balay } 965d5b485abSSatish Balay r_status3 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status)); 966d5b485abSSatish Balay CHKPTRQ(r_status3); 967d5b485abSSatish Balay s_status3 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 968d5b485abSSatish Balay CHKPTRQ(s_status3); 969d5b485abSSatish Balay 970d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 971d5b485abSSatish Balay sbuf_aa = (Scalar **)PetscMalloc((nrqr+1)*sizeof(Scalar *));CHKPTRQ(sbuf_aa); 972d5b485abSSatish Balay for ( i=0,j=0; i<nrqr; i++ ) j += req_size[i]; 973a2feddc5SSatish Balay sbuf_aa[0] = (Scalar*) PetscMalloc((j+1)*bs2*sizeof(Scalar));CHKPTRQ(sbuf_aa[0]); 974a2feddc5SSatish Balay for ( i=1; i<nrqr; i++ ) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 975d5b485abSSatish Balay 976d5b485abSSatish Balay s_waits4 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request)); 977d5b485abSSatish Balay CHKPTRQ(s_waits4); 978d5b485abSSatish Balay { 979d5b485abSSatish Balay for ( i=0; i<nrqr; i++ ) { 980d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 981d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 982d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 983d5b485abSSatish Balay ct2 = 0; 984d5b485abSSatish Balay for ( j=1,max1=rbuf1_i[0]; j<=max1; j++ ) { 985d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 986d5b485abSSatish Balay for ( k=0; k<kmax; k++,ct1++ ) { 987d5b485abSSatish Balay row = rbuf1_i[ct1] - cstart; 988d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 989d5b485abSSatish Balay ncols = nzA + nzB; 990d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 991a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 992d5b485abSSatish Balay 993d5b485abSSatish Balay /* load the column values for this row into vals*/ 994*5b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 995d5b485abSSatish Balay for ( l=0; l<nzB; l++ ) { 996a2feddc5SSatish Balay if ((bmap[cworkB[l]]) < cstart) 997a2feddc5SSatish Balay PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(Scalar)); 998d5b485abSSatish Balay else break; 999d5b485abSSatish Balay } 1000d5b485abSSatish Balay imark = l; 1001a2feddc5SSatish Balay for ( l=0; l<nzA; l++ ) 1002a2feddc5SSatish Balay PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(Scalar)); 1003a2feddc5SSatish Balay for ( l=imark; l<nzB; l++ ) 1004a2feddc5SSatish Balay PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(Scalar)); 1005d5b485abSSatish Balay ct2 += ncols; 1006d5b485abSSatish Balay } 1007d5b485abSSatish Balay } 1008a2feddc5SSatish Balay MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source[i],tag+3,comm,s_waits4+i); 1009d5b485abSSatish Balay } 1010d5b485abSSatish Balay } 1011d5b485abSSatish Balay r_status4 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status)); 1012d5b485abSSatish Balay CHKPTRQ(r_status4); 1013d5b485abSSatish Balay s_status4 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status)); 1014d5b485abSSatish Balay CHKPTRQ(s_status4); 1015d5b485abSSatish Balay PetscFree(rbuf1); 1016d5b485abSSatish Balay 1017d5b485abSSatish Balay /* Form the matrix */ 1018d5b485abSSatish Balay /* create col map */ 1019d5b485abSSatish Balay { 1020d5b485abSSatish Balay int *icol_i; 1021d5b485abSSatish Balay 1022a2feddc5SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->Nbs*sizeof(int); 1023d5b485abSSatish Balay cmap = (int **)PetscMalloc(len); CHKPTRQ(cmap); 1024d5b485abSSatish Balay cmap[0] = (int *)(cmap + ismax); 1025a2feddc5SSatish Balay PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int)); 1026a2feddc5SSatish Balay for ( i=1; i<ismax; i++ ) { cmap[i] = cmap[i-1] + c->Nbs; } 1027d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1028d5b485abSSatish Balay jmax = ncol[i]; 1029d5b485abSSatish Balay icol_i = icol[i]; 1030d5b485abSSatish Balay cmap_i = cmap[i]; 1031d5b485abSSatish Balay for ( j=0; j<jmax; j++ ) { 1032d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1033d5b485abSSatish Balay } 1034d5b485abSSatish Balay } 1035d5b485abSSatish Balay } 1036d5b485abSSatish Balay 1037d5b485abSSatish Balay 1038d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1039d5b485abSSatish Balay for ( i=0,j=0; i<ismax; i++ ) { j += nrow[i]; } 1040d5b485abSSatish Balay len = (1+ismax)*sizeof(int *) + j*sizeof(int); 1041d5b485abSSatish Balay lens = (int **)PetscMalloc(len); CHKPTRQ(lens); 1042d5b485abSSatish Balay lens[0] = (int *)(lens + ismax); 1043d5b485abSSatish Balay PetscMemzero(lens[0], j*sizeof(int)); 1044d5b485abSSatish Balay for ( i=1; i<ismax; i++ ) { lens[i] = lens[i-1] + nrow[i-1]; } 1045d5b485abSSatish Balay 1046d5b485abSSatish Balay /* Update lens from local data */ 1047d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1048d5b485abSSatish Balay jmax = nrow[i]; 1049d5b485abSSatish Balay cmap_i = cmap[i]; 1050d5b485abSSatish Balay irow_i = irow[i]; 1051d5b485abSSatish Balay lens_i = lens[i]; 1052d5b485abSSatish Balay for ( j=0; j<jmax; j++ ) { 1053d5b485abSSatish Balay row = irow_i[j]; 1054d5b485abSSatish Balay proc = rtable[row]; 1055d5b485abSSatish Balay if (proc == rank) { 105636f4e84dSSatish Balay /* Get indices from matA and then from matB */ 105736f4e84dSSatish Balay row = row - rstart; 105836f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 105936f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 106036f4e84dSSatish Balay for ( k=0; k<nzA; k++ ) 106136f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++;} 106236f4e84dSSatish Balay for ( k=0; k<nzB; k++ ) 106336f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++;} 1064d5b485abSSatish Balay } 1065d5b485abSSatish Balay } 1066a2feddc5SSatish Balay } 1067d5b485abSSatish Balay 1068d5b485abSSatish Balay /* Create row map*/ 1069a2feddc5SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->Mbs*sizeof(int); 1070d5b485abSSatish Balay rmap = (int **)PetscMalloc(len); CHKPTRQ(rmap); 1071d5b485abSSatish Balay rmap[0] = (int *)(rmap + ismax); 1072a2feddc5SSatish Balay PetscMemzero(rmap[0],ismax*c->Mbs*sizeof(int)); 1073a2feddc5SSatish Balay for ( i=1; i<ismax; i++ ) { rmap[i] = rmap[i-1] + c->Mbs;} 1074d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1075d5b485abSSatish Balay rmap_i = rmap[i]; 1076d5b485abSSatish Balay irow_i = irow[i]; 1077d5b485abSSatish Balay jmax = nrow[i]; 1078d5b485abSSatish Balay for ( j=0; j<jmax; j++ ) { 1079d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1080d5b485abSSatish Balay } 1081d5b485abSSatish Balay } 1082d5b485abSSatish Balay 1083d5b485abSSatish Balay /* Update lens from offproc data */ 1084d5b485abSSatish Balay { 1085d5b485abSSatish Balay int *rbuf2_i, *rbuf3_i, *sbuf1_i; 1086d5b485abSSatish Balay 1087d5b485abSSatish Balay for ( tmp2=0; tmp2<nrqs; tmp2++ ) { 1088d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits3, &i, r_status3+tmp2); 1089d5b485abSSatish Balay index = pa[i]; 1090d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1091d5b485abSSatish Balay jmax = sbuf1_i[0]; 1092d5b485abSSatish Balay ct1 = 2*jmax+1; 1093d5b485abSSatish Balay ct2 = 0; 1094d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1095d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1096d5b485abSSatish Balay for ( j=1; j<=jmax; j++ ) { 1097d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1098d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1099d5b485abSSatish Balay lens_i = lens[is_no]; 1100d5b485abSSatish Balay cmap_i = cmap[is_no]; 1101d5b485abSSatish Balay rmap_i = rmap[is_no]; 1102d5b485abSSatish Balay for ( k=0; k<max1; k++,ct1++ ) { 1103d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1104d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1105d5b485abSSatish Balay for ( l=0; l<max2; l++,ct2++ ) { 1106d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1107d5b485abSSatish Balay lens_i[row]++; 1108d5b485abSSatish Balay } 1109d5b485abSSatish Balay } 1110d5b485abSSatish Balay } 1111d5b485abSSatish Balay } 1112d5b485abSSatish Balay } 1113d5b485abSSatish Balay } 1114d5b485abSSatish Balay PetscFree(r_status3); PetscFree(r_waits3); 1115d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits3,s_status3); 1116d5b485abSSatish Balay PetscFree(s_status3); PetscFree(s_waits3); 1117d5b485abSSatish Balay 1118d5b485abSSatish Balay /* Create the submatrices */ 1119d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1120d5b485abSSatish Balay /* 1121d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1122d5b485abSSatish Balay */ 1123d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1124df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 112536f4e84dSSatish Balay if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) { 1126df36731bSBarry Smith SETERRQ(1,"MatGetSubmatrices_MPIBAIJ:Cannot reuse matrix. wrong size"); 1127d5b485abSSatish Balay } 112836f4e84dSSatish Balay if (PetscMemcmp(mat->ilen,lens[i], mat->mbs *sizeof(int))) { 1129df36731bSBarry Smith SETERRQ(1,"MatGetSubmatrices_MPIBAIJ:Cannot reuse matrix. wrong no of nonzeros"); 1130d5b485abSSatish Balay } 1131d5b485abSSatish Balay /* Initial matrix as if empty */ 113236f4e84dSSatish Balay PetscMemzero(mat->ilen,mat->mbs*sizeof(int)); 1133d5b485abSSatish Balay } 1134d5b485abSSatish Balay } 1135d5b485abSSatish Balay else { 1136d5b485abSSatish Balay *submat = submats = (Mat *)PetscMalloc(ismax*sizeof(Mat)); CHKPTRQ(submats); 1137d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1138a2feddc5SSatish Balay ierr = MatCreateSeqBAIJ(comm,a->bs,nrow[i]*bs,ncol[i]*bs,0,lens[i],submats+i);CHKERRQ(ierr); 1139d5b485abSSatish Balay } 1140d5b485abSSatish Balay } 1141d5b485abSSatish Balay 1142d5b485abSSatish Balay /* Assemble the matrices */ 1143d5b485abSSatish Balay /* First assemble the local rows */ 1144d5b485abSSatish Balay { 114536f4e84dSSatish Balay int ilen_row,*imat_ilen, *imat_j, *imat_i; 1146d5b485abSSatish Balay Scalar *imat_a; 1147d5b485abSSatish Balay 1148d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1149df36731bSBarry Smith mat = (Mat_SeqBAIJ *) submats[i]->data; 1150d5b485abSSatish Balay imat_ilen = mat->ilen; 1151d5b485abSSatish Balay imat_j = mat->j; 1152d5b485abSSatish Balay imat_i = mat->i; 1153d5b485abSSatish Balay imat_a = mat->a; 1154d5b485abSSatish Balay cmap_i = cmap[i]; 1155d5b485abSSatish Balay rmap_i = rmap[i]; 1156d5b485abSSatish Balay irow_i = irow[i]; 1157d5b485abSSatish Balay jmax = nrow[i]; 1158d5b485abSSatish Balay for ( j=0; j<jmax; j++ ) { 1159d5b485abSSatish Balay row = irow_i[j]; 1160d5b485abSSatish Balay proc = rtable[row]; 1161d5b485abSSatish Balay if (proc == rank) { 116236f4e84dSSatish Balay row = row - rstart; 116336f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 116436f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 116536f4e84dSSatish Balay cworkA = a_j + a_i[row]; 116636f4e84dSSatish Balay cworkB = b_j + b_i[row]; 116736f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 116836f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 116936f4e84dSSatish Balay 117036f4e84dSSatish Balay row = rmap_i[row + rstart]; 1171df36731bSBarry Smith mat_i = imat_i[row]; 117236f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1173d5b485abSSatish Balay mat_j = imat_j + mat_i; 117436f4e84dSSatish Balay ilen_row = imat_ilen[row]; 117536f4e84dSSatish Balay 117636f4e84dSSatish Balay /* load the column indices for this row into cols*/ 117736f4e84dSSatish Balay for ( l=0; l<nzB; l++ ) { 117836f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 117936f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1180df36731bSBarry Smith *mat_j++ = tcol - 1; 118136f4e84dSSatish Balay PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 1182d5b485abSSatish Balay ilen_row++; 1183d5b485abSSatish Balay } 1184d5b485abSSatish Balay } 118536f4e84dSSatish Balay else break; 118636f4e84dSSatish Balay } 118736f4e84dSSatish Balay imark = l; 118836f4e84dSSatish Balay for ( l=0; l<nzA; l++ ) { 118936f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 119036f4e84dSSatish Balay *mat_j++ = tcol - 1; 119136f4e84dSSatish Balay PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 119236f4e84dSSatish Balay ilen_row++; 119336f4e84dSSatish Balay } 119436f4e84dSSatish Balay } 119536f4e84dSSatish Balay for ( l=imark; l<nzB; l++ ) { 119636f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 119736f4e84dSSatish Balay *mat_j++ = tcol - 1; 119836f4e84dSSatish Balay PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 119936f4e84dSSatish Balay ilen_row++; 120036f4e84dSSatish Balay } 120136f4e84dSSatish Balay } 1202d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1203d5b485abSSatish Balay } 1204d5b485abSSatish Balay } 120536f4e84dSSatish Balay 1206d5b485abSSatish Balay } 1207d5b485abSSatish Balay } 1208d5b485abSSatish Balay 1209d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1210d5b485abSSatish Balay { 1211d5b485abSSatish Balay int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1212d5b485abSSatish Balay int *imat_j,*imat_i; 1213d5b485abSSatish Balay Scalar *imat_a,*rbuf4_i; 1214d5b485abSSatish Balay 1215d5b485abSSatish Balay for ( tmp2=0; tmp2<nrqs; tmp2++ ) { 1216d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits4, &i, r_status4+tmp2); 1217d5b485abSSatish Balay index = pa[i]; 1218d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1219d5b485abSSatish Balay jmax = sbuf1_i[0]; 1220d5b485abSSatish Balay ct1 = 2*jmax + 1; 1221d5b485abSSatish Balay ct2 = 0; 1222d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1223d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1224d5b485abSSatish Balay rbuf4_i = rbuf4[i]; 1225d5b485abSSatish Balay for ( j=1; j<=jmax; j++ ) { 1226d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1227d5b485abSSatish Balay rmap_i = rmap[is_no]; 1228d5b485abSSatish Balay cmap_i = cmap[is_no]; 1229df36731bSBarry Smith mat = (Mat_SeqBAIJ *) submats[is_no]->data; 1230d5b485abSSatish Balay imat_ilen = mat->ilen; 1231d5b485abSSatish Balay imat_j = mat->j; 1232d5b485abSSatish Balay imat_i = mat->i; 1233d5b485abSSatish Balay imat_a = mat->a; 1234d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1235d5b485abSSatish Balay for ( k=0; k<max1; k++, ct1++ ) { 1236d5b485abSSatish Balay row = sbuf1_i[ct1]; 1237d5b485abSSatish Balay row = rmap_i[row]; 1238d5b485abSSatish Balay ilen = imat_ilen[row]; 1239df36731bSBarry Smith mat_i = imat_i[row]; 124036f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1241d5b485abSSatish Balay mat_j = imat_j + mat_i; 1242d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1243d5b485abSSatish Balay for ( l=0; l<max2; l++,ct2++ ) { 1244d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1245df36731bSBarry Smith *mat_j++ = tcol - 1; 124636f4e84dSSatish Balay /* *mat_a++ = rbuf4_i[ct2]; */ 124736f4e84dSSatish Balay PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 1248d5b485abSSatish Balay ilen++; 1249d5b485abSSatish Balay } 1250d5b485abSSatish Balay } 1251d5b485abSSatish Balay imat_ilen[row] = ilen; 1252d5b485abSSatish Balay } 1253d5b485abSSatish Balay } 1254d5b485abSSatish Balay } 1255d5b485abSSatish Balay } 1256d5b485abSSatish Balay PetscFree(r_status4); PetscFree(r_waits4); 1257d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits4,s_status4); 1258d5b485abSSatish Balay PetscFree(s_waits4); PetscFree(s_status4); 1259d5b485abSSatish Balay 1260d5b485abSSatish Balay /* Restore the indices */ 1261d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 1262d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i], irow+i); CHKERRQ(ierr); 1263d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i], icol+i); CHKERRQ(ierr); 1264d5b485abSSatish Balay } 1265d5b485abSSatish Balay 1266d5b485abSSatish Balay /* Destroy allocated memory */ 1267d5b485abSSatish Balay PetscFree(irow); 1268d5b485abSSatish Balay PetscFree(w1); 1269d5b485abSSatish Balay PetscFree(pa); 1270d5b485abSSatish Balay 1271d5b485abSSatish Balay PetscFree(sbuf1); 1272d5b485abSSatish Balay PetscFree(rbuf2); 1273d5b485abSSatish Balay for ( i=0; i<nrqr; ++i ) { 1274d5b485abSSatish Balay PetscFree(sbuf2[i]); 1275d5b485abSSatish Balay } 1276d5b485abSSatish Balay for ( i=0; i<nrqs; ++i ) { 1277d5b485abSSatish Balay PetscFree(rbuf3[i]); 1278d5b485abSSatish Balay PetscFree(rbuf4[i]); 1279d5b485abSSatish Balay } 1280d5b485abSSatish Balay 1281d5b485abSSatish Balay PetscFree(sbuf2); 1282d5b485abSSatish Balay PetscFree(rbuf3); 1283d5b485abSSatish Balay PetscFree(rbuf4 ); 1284d5b485abSSatish Balay PetscFree(sbuf_aj[0]); 1285d5b485abSSatish Balay PetscFree(sbuf_aj); 1286d5b485abSSatish Balay PetscFree(sbuf_aa[0]); 1287d5b485abSSatish Balay PetscFree(sbuf_aa); 1288d5b485abSSatish Balay 1289d5b485abSSatish Balay PetscFree(cmap); 1290d5b485abSSatish Balay PetscFree(rmap); 1291d5b485abSSatish Balay PetscFree(lens); 1292d5b485abSSatish Balay 1293d5b485abSSatish Balay for ( i=0; i<ismax; i++ ) { 129436f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 129536f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1296d5b485abSSatish Balay } 1297d5b485abSSatish Balay return 0; 1298d5b485abSSatish Balay } 1299