1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3a40ed3dSBarry Smith static char vcid[] = "$Id: baijov.c,v 1.25 1997/10/06 16:24:17 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 */ 870f55243SBarry 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 185615d1e5SSatish Balay #undef __FUNC__ 19d4bb536fSBarry Smith #define __FUNC__ "MatCompressIndicesGeneral_MPIBAIJ" 20e757655aSSatish Balay static int MatCompressIndicesGeneral_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 21df36731bSBarry Smith { 22df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 23df36731bSBarry Smith int ierr,isz,bs = baij->bs,Nbs,n,i,j,*idx,*nidx,ival; 24bbd702dbSSatish Balay BT table; 25df36731bSBarry Smith 26*3a40ed3dSBarry Smith PetscFunctionBegin; 27df36731bSBarry Smith Nbs = baij->Nbs; 28df36731bSBarry Smith nidx = (int *) PetscMalloc((Nbs+1)*sizeof(int)); CHKPTRQ(nidx); 29bbd702dbSSatish Balay ierr = BTCreate(Nbs,table); CHKERRQ(ierr); 30df36731bSBarry Smith 31df36731bSBarry Smith for (i=0; i<imax; i++) { 32df36731bSBarry Smith isz = 0; 33bbd702dbSSatish Balay BTMemzero(Nbs,table); 3436f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 3536f4e84dSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 36e757655aSSatish Balay for (j=0; j<n ; j++) { 37df36731bSBarry Smith ival = idx[j]/bs; /* convert the indices into block indices */ 38e3372554SBarry Smith if (ival>Nbs) SETERRQ(1,0,"index greater than mat-dim"); 39bbd702dbSSatish Balay if(!BTLookupSet(table, ival)) { nidx[isz++] = ival;} 40df36731bSBarry Smith } 4136f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 42029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is_out+i)); CHKERRQ(ierr); 43df36731bSBarry Smith } 44bbd702dbSSatish Balay BTDestroy(table); 45df36731bSBarry Smith PetscFree(nidx); 46*3a40ed3dSBarry Smith PetscFunctionReturn(0); 47df36731bSBarry Smith } 48df36731bSBarry Smith 495615d1e5SSatish Balay #undef __FUNC__ 50d4bb536fSBarry Smith #define __FUNC__ "MatCompressIndicesSorted_MPIBAIJ" 51e757655aSSatish Balay static int MatCompressIndicesSorted_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 52e757655aSSatish Balay { 53e757655aSSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 54e757655aSSatish Balay int ierr,bs=baij->bs,i,j,k,val,n,*idx,*nidx,Nbs=baij->Nbs,*idx_local; 55e757655aSSatish Balay PetscTruth flg; 56e757655aSSatish Balay 57*3a40ed3dSBarry Smith PetscFunctionBegin; 58e757655aSSatish Balay for (i=0; i<imax; i++) { 59e757655aSSatish Balay ierr = ISSorted(is_in[i],&flg); CHKERRQ(ierr); 60e3372554SBarry Smith if (!flg) SETERRQ(1,0,"Indices are not sorted"); 61e757655aSSatish Balay } 62e757655aSSatish Balay nidx = (int *) PetscMalloc((Nbs+1)*sizeof(int)); CHKPTRQ(nidx); 63*3a40ed3dSBarry Smith /* Now check if the indices are in block order */ 64e757655aSSatish Balay for (i=0; i<imax; i++) { 65e757655aSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 66e757655aSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 67e3372554SBarry Smith if ( n%bs !=0 ) SETERRA(1,0,"Indices are not block ordered"); 68e757655aSSatish Balay 69e757655aSSatish Balay n = n/bs; /* The reduced index size */ 70e757655aSSatish Balay idx_local = idx; 71e757655aSSatish Balay for (j=0; j<n ; j++) { 72e757655aSSatish Balay val = idx_local[0]; 73e3372554SBarry Smith if(val%bs != 0) SETERRA(1,0,"Indices are not block ordered"); 74e757655aSSatish Balay for (k=0; k<bs; k++) { 75e3372554SBarry Smith if ( val+k != idx_local[k]) SETERRA(1,0,"Indices are not block ordered"); 76e757655aSSatish Balay } 77e757655aSSatish Balay nidx[j] = val/bs; 78e757655aSSatish Balay idx_local +=bs; 79e757655aSSatish Balay } 80e757655aSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 81029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,nidx,(is_out+i)); CHKERRQ(ierr); 82e757655aSSatish Balay } 83e757655aSSatish Balay PetscFree(nidx); 84*3a40ed3dSBarry Smith PetscFunctionReturn(0); 85e757655aSSatish Balay } 86e757655aSSatish Balay 875615d1e5SSatish Balay #undef __FUNC__ 88d4bb536fSBarry Smith #define __FUNC__ "MatExpandIndices_MPIBAIJ" 8936f4e84dSSatish Balay static int MatExpandIndices_MPIBAIJ(Mat C, int imax, IS *is_in, IS *is_out) 90df36731bSBarry Smith { 91df36731bSBarry Smith Mat_MPIBAIJ *baij = (Mat_MPIBAIJ *) C->data; 92df36731bSBarry Smith int ierr,bs = baij->bs,Nbs,n,i,j,k,*idx,*nidx; 93df36731bSBarry Smith 94*3a40ed3dSBarry Smith PetscFunctionBegin; 95df36731bSBarry Smith Nbs = baij->Nbs; 96df36731bSBarry Smith 97df36731bSBarry Smith nidx = (int *) PetscMalloc((Nbs*bs+1)*sizeof(int)); CHKPTRQ(nidx); 98df36731bSBarry Smith 99df36731bSBarry Smith for ( i=0; i<imax; i++ ) { 10036f4e84dSSatish Balay ierr = ISGetIndices(is_in[i],&idx); CHKERRQ(ierr); 10136f4e84dSSatish Balay ierr = ISGetSize(is_in[i],&n); CHKERRQ(ierr); 102df36731bSBarry Smith for (j=0; j<n ; ++j){ 103df36731bSBarry Smith for (k=0; k<bs; k++) 104df36731bSBarry Smith nidx[j*bs+k] = idx[j]*bs+k; 105df36731bSBarry Smith } 10636f4e84dSSatish Balay ierr = ISRestoreIndices(is_in[i],&idx); CHKERRQ(ierr); 107029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, n*bs, nidx, (is_out+i)); CHKERRQ(ierr); 108df36731bSBarry Smith } 109df36731bSBarry Smith PetscFree(nidx); 110*3a40ed3dSBarry Smith PetscFunctionReturn(0); 111df36731bSBarry Smith } 112df36731bSBarry Smith 113df36731bSBarry Smith 1145615d1e5SSatish Balay #undef __FUNC__ 115d4bb536fSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ" 116df36731bSBarry Smith int MatIncreaseOverlap_MPIBAIJ(Mat C, int imax, IS *is, int ov) 117d5b485abSSatish Balay { 118d5b485abSSatish Balay int i, ierr; 11936f4e84dSSatish Balay IS *is_new; 12036f4e84dSSatish Balay 121*3a40ed3dSBarry Smith PetscFunctionBegin; 12236f4e84dSSatish Balay is_new = (IS *)PetscMalloc(imax*sizeof(IS)); CHKPTRQ(is_new); 123df36731bSBarry Smith /* Convert the indices into block format */ 124e757655aSSatish Balay ierr = MatCompressIndicesGeneral_MPIBAIJ(C, imax, is,is_new); CHKERRQ(ierr); 125e3372554SBarry Smith if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative overlap specified\n");} 126d5b485abSSatish Balay for (i=0; i<ov; ++i) { 12736f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new); CHKERRQ(ierr); 128d5b485abSSatish Balay } 12936f4e84dSSatish Balay for (i=0; i<imax; i++) ISDestroy(is[i]); 13036f4e84dSSatish Balay ierr = MatExpandIndices_MPIBAIJ(C, imax, is_new,is); CHKERRQ(ierr); 13136f4e84dSSatish Balay for (i=0; i<imax; i++) ISDestroy(is_new[i]); 13236f4e84dSSatish Balay PetscFree(is_new); 133*3a40ed3dSBarry Smith PetscFunctionReturn(0); 134d5b485abSSatish Balay } 135d5b485abSSatish Balay 136d5b485abSSatish Balay /* 137d5b485abSSatish Balay Sample message format: 138d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 139d5b485abSSatish Balay to index sets 1s[1], is[5] 140d5b485abSSatish Balay mesg [0] = 2 ( no of index sets in the mesg) 141d5b485abSSatish Balay ----------- 142d5b485abSSatish Balay mesg [1] = 1 => is[1] 143d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 144d5b485abSSatish Balay ----------- 145d5b485abSSatish Balay mesg [5] = 5 => is[5] 146d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 147d5b485abSSatish Balay ----------- 148d5b485abSSatish Balay mesg [7] 149d5b485abSSatish Balay mesg [n] datas[1] 150d5b485abSSatish Balay ----------- 151d5b485abSSatish Balay mesg[n+1] 152d5b485abSSatish Balay mesg[m] data(is[5]) 153d5b485abSSatish Balay ----------- 154d5b485abSSatish Balay 155d5b485abSSatish Balay Notes: 156d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 157d5b485abSSatish Balay nrqr - no of requests recieved (which have to be or which have been processed 158d5b485abSSatish Balay */ 1595615d1e5SSatish Balay #undef __FUNC__ 160d4bb536fSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Once" 161df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Once(Mat C, int imax, IS *is) 162d5b485abSSatish Balay { 163df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 164d5b485abSSatish Balay int **idx, *n, *w1, *w2, *w3, *w4, *rtable,**data,len,*idx_i; 165df36731bSBarry Smith int size,rank,Mbs,i,j,k,ierr,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 166d5b485abSSatish Balay int *ctr,*pa,tag,*tmp,bsz,nrqr,*isz,*isz1,**xdata,bsz1,**rbuf2; 167bbd702dbSSatish Balay BT *table; 168d5b485abSSatish Balay MPI_Comm comm; 169d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 170d5b485abSSatish Balay MPI_Status *s_status,*recv_status; 171d5b485abSSatish Balay 172*3a40ed3dSBarry Smith PetscFunctionBegin; 173d5b485abSSatish Balay comm = C->comm; 174d5b485abSSatish Balay tag = C->tag; 175d5b485abSSatish Balay size = c->size; 176d5b485abSSatish Balay rank = c->rank; 177df36731bSBarry Smith Mbs = c->Mbs; 178d5b485abSSatish Balay 179df36731bSBarry Smith len = (imax+1)*sizeof(int *) + (imax + Mbs)*sizeof(int); 180d5b485abSSatish Balay idx = (int **) PetscMalloc(len); CHKPTRQ(idx); 181d5b485abSSatish Balay n = (int *) (idx + imax); 182d5b485abSSatish Balay rtable = n + imax; 183d5b485abSSatish Balay 184d5b485abSSatish Balay for (i=0; i<imax; i++) { 185d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]); CHKERRQ(ierr); 186d5b485abSSatish Balay ierr = ISGetSize(is[i],&n[i]); CHKERRQ(ierr); 187d5b485abSSatish Balay } 188d5b485abSSatish Balay 189d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 190d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 191d5b485abSSatish Balay len = c->rowners[i+1]; 192d5b485abSSatish Balay for (; j<len; j++) { 193d5b485abSSatish Balay rtable[j] = i; 194d5b485abSSatish Balay } 195d5b485abSSatish Balay } 196d5b485abSSatish Balay 197d5b485abSSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 198d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 199d5b485abSSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1);/* mesg size */ 200d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 201d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 202d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 203d5b485abSSatish Balay PetscMemzero(w1,size*3*sizeof(int)); /* initialise work vector*/ 204d5b485abSSatish Balay for (i=0; i<imax; i++) { 205d5b485abSSatish Balay PetscMemzero(w4,size*sizeof(int)); /* initialise work vector*/ 206d5b485abSSatish Balay idx_i = idx[i]; 207d5b485abSSatish Balay len = n[i]; 208d5b485abSSatish Balay for (j=0; j<len; j++) { 209d5b485abSSatish Balay row = idx_i[j]; 210d5b485abSSatish Balay proc = rtable[row]; 211d5b485abSSatish Balay w4[proc]++; 212d5b485abSSatish Balay } 213d5b485abSSatish Balay for (j=0; j<size; j++){ 214d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 215d5b485abSSatish Balay } 216d5b485abSSatish Balay } 217d5b485abSSatish Balay 218d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 219d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 220d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 221d5b485abSSatish Balay w3[rank] = 0; 222d5b485abSSatish Balay for ( i=0; i<size; i++) { 223d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 224d5b485abSSatish Balay } 225d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 226d5b485abSSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); 227d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 228d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 229d5b485abSSatish Balay } 230d5b485abSSatish Balay 231d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 232d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 233d5b485abSSatish Balay j = pa[i]; 234d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 235d5b485abSSatish Balay msz += w1[j]; 236d5b485abSSatish Balay } 237d5b485abSSatish Balay 238d5b485abSSatish Balay 239d5b485abSSatish Balay /* Do a global reduction to determine how many messages to expect*/ 240d5b485abSSatish Balay { 241d5b485abSSatish Balay int *rw1, *rw2; 242d5b485abSSatish Balay rw1 = (int *) PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 243d5b485abSSatish Balay rw2 = rw1+size; 244d5b485abSSatish Balay MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm); 245d5b485abSSatish Balay bsz = rw1[rank]; 246d5b485abSSatish Balay MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm); 247d5b485abSSatish Balay nrqr = rw2[rank]; 248d5b485abSSatish Balay PetscFree(rw1); 249d5b485abSSatish Balay } 250d5b485abSSatish Balay 251d5b485abSSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 252d5b485abSSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 253d5b485abSSatish Balay rbuf = (int**) PetscMalloc(len); CHKPTRQ(rbuf); 254d5b485abSSatish Balay rbuf[0] = (int *) (rbuf + nrqr); 255d5b485abSSatish Balay for (i=1; i<nrqr; ++i) rbuf[i] = rbuf[i-1] + bsz; 256d5b485abSSatish Balay 257d5b485abSSatish Balay /* Post the receives */ 258*3a40ed3dSBarry Smith r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1); 259d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 260d5b485abSSatish Balay MPI_Irecv(rbuf[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i); 261d5b485abSSatish Balay } 262d5b485abSSatish Balay 263d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 264d5b485abSSatish Balay len = 2*size*sizeof(int*) + (size+msz)*sizeof(int); 265d5b485abSSatish Balay outdat = (int **)PetscMalloc(len); CHKPTRQ(outdat); 266d5b485abSSatish Balay ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 267d5b485abSSatish Balay PetscMemzero(outdat,2*size*sizeof(int*)); 268d5b485abSSatish Balay tmp = (int *) (outdat + 2*size); 269d5b485abSSatish Balay ctr = tmp + msz; 270d5b485abSSatish Balay 271d5b485abSSatish Balay { 272d5b485abSSatish Balay int *iptr = tmp,ict = 0; 273d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 274d5b485abSSatish Balay j = pa[i]; 275d5b485abSSatish Balay iptr += ict; 276d5b485abSSatish Balay outdat[j] = iptr; 277d5b485abSSatish Balay ict = w1[j]; 278d5b485abSSatish Balay } 279d5b485abSSatish Balay } 280d5b485abSSatish Balay 281d5b485abSSatish Balay /* Form the outgoing messages */ 282d5b485abSSatish Balay /*plug in the headers*/ 283d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 284d5b485abSSatish Balay j = pa[i]; 285d5b485abSSatish Balay outdat[j][0] = 0; 286d5b485abSSatish Balay PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int)); 287d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 288d5b485abSSatish Balay } 289d5b485abSSatish Balay 290d5b485abSSatish Balay /* Memory for doing local proc's work*/ 291d5b485abSSatish Balay { 292d5b485abSSatish Balay int *d_p; 293d5b485abSSatish Balay char *t_p; 294d5b485abSSatish Balay 295bbd702dbSSatish Balay len = (imax)*(sizeof(BT) + sizeof(int *) + sizeof(int)) + 2962a8f2162SSatish Balay (Mbs)*imax*sizeof(int) + (Mbs/BITSPERBYTE+1)*imax*sizeof(char) + 1; 297bbd702dbSSatish Balay table = (BT *)PetscMalloc(len); CHKPTRQ(table); 298bbd702dbSSatish Balay PetscMemzero(table,len); 299d5b485abSSatish Balay data = (int **)(table + imax); 300bbd702dbSSatish Balay isz = (int *)(data + imax); 301bbd702dbSSatish Balay d_p = (int *)(isz + imax); 302bbd702dbSSatish Balay t_p = (char *)(d_p + Mbs*imax); 303bbd702dbSSatish Balay for (i=0; i<imax; i++) { 3042a8f2162SSatish Balay table[i] = t_p + (Mbs/BITSPERBYTE+1)*i; 305bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 306d5b485abSSatish Balay } 307d5b485abSSatish Balay } 308d5b485abSSatish Balay 309d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 310d5b485abSSatish Balay { 311d5b485abSSatish Balay int n_i,*data_i,isz_i,*outdat_j,ctr_j; 312bbd702dbSSatish Balay BT table_i; 313d5b485abSSatish Balay 314d5b485abSSatish Balay for (i=0; i<imax; i++) { 315d5b485abSSatish Balay PetscMemzero(ctr,size*sizeof(int)); 316d5b485abSSatish Balay n_i = n[i]; 317d5b485abSSatish Balay table_i = table[i]; 318d5b485abSSatish Balay idx_i = idx[i]; 319d5b485abSSatish Balay data_i = data[i]; 320d5b485abSSatish Balay isz_i = isz[i]; 321d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 322d5b485abSSatish Balay row = idx_i[j]; 323d5b485abSSatish Balay proc = rtable[row]; 324d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 325d5b485abSSatish Balay ctr[proc]++; 326d5b485abSSatish Balay *ptr[proc] = row; 327d5b485abSSatish Balay ptr[proc]++; 328d5b485abSSatish Balay } 329d5b485abSSatish Balay else { /* Update the local table */ 330bbd702dbSSatish Balay if (!BTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 331d5b485abSSatish Balay } 332d5b485abSSatish Balay } 333d5b485abSSatish Balay /* Update the headers for the current IS */ 334d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 335d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 336d5b485abSSatish Balay outdat_j = outdat[j]; 337d5b485abSSatish Balay k = ++outdat_j[0]; 338d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 339d5b485abSSatish Balay outdat_j[2*k-1] = i; 340d5b485abSSatish Balay } 341d5b485abSSatish Balay } 342d5b485abSSatish Balay isz[i] = isz_i; 343d5b485abSSatish Balay } 344d5b485abSSatish Balay } 345d5b485abSSatish Balay 346d5b485abSSatish Balay 347d5b485abSSatish Balay 348d5b485abSSatish Balay /* Now post the sends */ 349*3a40ed3dSBarry Smith s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1); 350d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 351d5b485abSSatish Balay j = pa[i]; 352d5b485abSSatish Balay MPI_Isend(outdat[j], w1[j], MPI_INT, j, tag, comm, s_waits1+i); 353d5b485abSSatish Balay } 354d5b485abSSatish Balay 355d5b485abSSatish Balay /* No longer need the original indices*/ 356d5b485abSSatish Balay for (i=0; i<imax; ++i) { 357d5b485abSSatish Balay ierr = ISRestoreIndices(is[i], idx+i); CHKERRQ(ierr); 358d5b485abSSatish Balay } 359d5b485abSSatish Balay PetscFree(idx); 360d5b485abSSatish Balay 361d5b485abSSatish Balay for (i=0; i<imax; ++i) { 362d5b485abSSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 363d5b485abSSatish Balay } 364d5b485abSSatish Balay 365d5b485abSSatish Balay /* Do Local work*/ 366df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 367d5b485abSSatish Balay 368d5b485abSSatish Balay /* Receive messages*/ 369d5b485abSSatish Balay { 370d5b485abSSatish Balay int index; 371d5b485abSSatish Balay 372*3a40ed3dSBarry Smith recv_status = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(recv_status); 373d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 374d5b485abSSatish Balay MPI_Waitany(nrqr, r_waits1, &index, recv_status+i); 375d5b485abSSatish Balay } 376d5b485abSSatish Balay 377*3a40ed3dSBarry Smith s_status = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status); 378d5b485abSSatish Balay MPI_Waitall(nrqs,s_waits1,s_status); 379d5b485abSSatish Balay } 380d5b485abSSatish Balay 381d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 382d5b485abSSatish Balay PetscFree(outdat); 383d5b485abSSatish Balay PetscFree(w1); 384d5b485abSSatish Balay 385d5b485abSSatish Balay xdata = (int **)PetscMalloc((nrqr+1)*sizeof(int *)); CHKPTRQ(xdata); 386d5b485abSSatish Balay isz1 = (int *)PetscMalloc((nrqr+1)*sizeof(int)); CHKPTRQ(isz1); 387df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 388d5b485abSSatish Balay PetscFree(rbuf); 389d5b485abSSatish Balay 390d5b485abSSatish Balay /* Send the data back*/ 391d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 392d5b485abSSatish Balay { 393d5b485abSSatish Balay int *rw1, *rw2; 394d5b485abSSatish Balay 395d5b485abSSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 396d5b485abSSatish Balay PetscMemzero(rw1,2*size*sizeof(int)); 397d5b485abSSatish Balay rw2 = rw1+size; 398d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 399d5b485abSSatish Balay proc = recv_status[i].MPI_SOURCE; 400d5b485abSSatish Balay rw1[proc] = isz1[i]; 401d5b485abSSatish Balay } 402d5b485abSSatish Balay 403d5b485abSSatish Balay MPI_Allreduce(rw1, rw2, size, MPI_INT, MPI_MAX, comm); 404d5b485abSSatish Balay bsz1 = rw2[rank]; 405d5b485abSSatish Balay PetscFree(rw1); 406d5b485abSSatish Balay } 407d5b485abSSatish Balay 408d5b485abSSatish Balay /* Allocate buffers*/ 409d5b485abSSatish Balay 410d5b485abSSatish Balay /* Allocate memory for recv buffers. Prob none if nrqr = 0 ???? */ 411d5b485abSSatish Balay len = (nrqs+1)*sizeof(int*) + nrqs*bsz1*sizeof(int); 412d5b485abSSatish Balay rbuf2 = (int**) PetscMalloc(len); CHKPTRQ(rbuf2); 413d5b485abSSatish Balay rbuf2[0] = (int *) (rbuf2 + nrqs); 414d5b485abSSatish Balay for (i=1; i<nrqs; ++i) rbuf2[i] = rbuf2[i-1] + bsz1; 415d5b485abSSatish Balay 416d5b485abSSatish Balay /* Post the receives */ 417*3a40ed3dSBarry Smith r_waits2 = (MPI_Request *)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2); 418d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 419d5b485abSSatish Balay MPI_Irecv(rbuf2[i], bsz1, MPI_INT, MPI_ANY_SOURCE, tag, comm, r_waits2+i); 420d5b485abSSatish Balay } 421d5b485abSSatish Balay 422d5b485abSSatish Balay /* Now post the sends */ 423*3a40ed3dSBarry Smith s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2); 424d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 425d5b485abSSatish Balay j = recv_status[i].MPI_SOURCE; 426d5b485abSSatish Balay MPI_Isend( xdata[i], isz1[i], MPI_INT, j, tag, comm, s_waits2+i); 427d5b485abSSatish Balay } 428d5b485abSSatish Balay 429d5b485abSSatish Balay /* receive work done on other processors*/ 430d5b485abSSatish Balay { 431d5b485abSSatish Balay int index, is_no, ct1, max,*rbuf2_i,isz_i,*data_i,jmax; 432bbd702dbSSatish Balay BT table_i; 433d5b485abSSatish Balay MPI_Status *status2; 434d5b485abSSatish Balay 435d5b485abSSatish Balay status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(status2); 436d5b485abSSatish Balay 437d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 438d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits2, &index, status2+i); 439d5b485abSSatish Balay /* Process the message*/ 440d5b485abSSatish Balay rbuf2_i = rbuf2[index]; 441d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 442d5b485abSSatish Balay jmax = rbuf2[index][0]; 443d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 444d5b485abSSatish Balay max = rbuf2_i[2*j]; 445d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 446d5b485abSSatish Balay isz_i = isz[is_no]; 447d5b485abSSatish Balay data_i = data[is_no]; 448d5b485abSSatish Balay table_i = table[is_no]; 449d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 450d5b485abSSatish Balay row = rbuf2_i[ct1]; 451bbd702dbSSatish Balay if (!BTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 452d5b485abSSatish Balay } 453d5b485abSSatish Balay isz[is_no] = isz_i; 454d5b485abSSatish Balay } 455d5b485abSSatish Balay } 456d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits2,status2); 457d5b485abSSatish Balay PetscFree(status2); 458d5b485abSSatish Balay } 459d5b485abSSatish Balay 460d5b485abSSatish Balay for (i=0; i<imax; ++i) { 461029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz[i], data[i], is+i); CHKERRQ(ierr); 462d5b485abSSatish Balay } 463d5b485abSSatish Balay 464d5b485abSSatish Balay PetscFree(pa); 465d5b485abSSatish Balay PetscFree(rbuf2); 466d5b485abSSatish Balay PetscFree(s_waits1); 467d5b485abSSatish Balay PetscFree(r_waits1); 468d5b485abSSatish Balay PetscFree(s_waits2); 469d5b485abSSatish Balay PetscFree(r_waits2); 470d5b485abSSatish Balay PetscFree(table); 471d5b485abSSatish Balay PetscFree(s_status); 472d5b485abSSatish Balay PetscFree(recv_status); 473d5b485abSSatish Balay PetscFree(xdata[0]); 474d5b485abSSatish Balay PetscFree(xdata); 475d5b485abSSatish Balay PetscFree(isz1); 476*3a40ed3dSBarry Smith PetscFunctionReturn(0); 477d5b485abSSatish Balay } 478d5b485abSSatish Balay 4795615d1e5SSatish Balay #undef __FUNC__ 480d4bb536fSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Local" 481d5b485abSSatish Balay /* 482df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 483d5b485abSSatish Balay the work on the local processor. 484d5b485abSSatish Balay 485d5b485abSSatish Balay Inputs: 486df36731bSBarry Smith C - MAT_MPIBAIJ; 487d5b485abSSatish Balay imax - total no of index sets processed at a time; 488df36731bSBarry Smith table - an array of char - size = Mbs bits. 489d5b485abSSatish Balay 490d5b485abSSatish Balay Output: 491d5b485abSSatish Balay isz - array containing the count of the solution elements correspondign 492d5b485abSSatish Balay to each index set; 493d5b485abSSatish Balay data - pointer to the solutions 494d5b485abSSatish Balay */ 495bbd702dbSSatish Balay static int MatIncreaseOverlap_MPIBAIJ_Local(Mat C,int imax,BT *table,int *isz, 496d5b485abSSatish Balay int **data) 497d5b485abSSatish Balay { 498df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 499d5b485abSSatish Balay Mat A = c->A, B = c->B; 500df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 501df36731bSBarry Smith int start, end, val, max, rstart,cstart,*ai, *aj; 502d5b485abSSatish Balay int *bi, *bj, *garray, i, j, k, row,*data_i,isz_i; 503bbd702dbSSatish Balay BT table_i; 504d5b485abSSatish Balay 505*3a40ed3dSBarry Smith PetscFunctionBegin; 506d5b485abSSatish Balay rstart = c->rstart; 507d5b485abSSatish Balay cstart = c->cstart; 508d5b485abSSatish Balay ai = a->i; 509df36731bSBarry Smith aj = a->j; 510d5b485abSSatish Balay bi = b->i; 511df36731bSBarry Smith bj = b->j; 512d5b485abSSatish Balay garray = c->garray; 513d5b485abSSatish Balay 514d5b485abSSatish Balay 515d5b485abSSatish Balay for (i=0; i<imax; i++) { 516d5b485abSSatish Balay data_i = data[i]; 517d5b485abSSatish Balay table_i = table[i]; 518d5b485abSSatish Balay isz_i = isz[i]; 519d5b485abSSatish Balay for (j=0, max=isz[i]; j<max; j++) { 520d5b485abSSatish Balay row = data_i[j] - rstart; 521d5b485abSSatish Balay start = ai[row]; 522d5b485abSSatish Balay end = ai[row+1]; 523d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 524df36731bSBarry Smith val = aj[k] + cstart; 525bbd702dbSSatish Balay if (!BTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 526d5b485abSSatish Balay } 527d5b485abSSatish Balay start = bi[row]; 528d5b485abSSatish Balay end = bi[row+1]; 529d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 530df36731bSBarry Smith val = garray[bj[k]]; 531bbd702dbSSatish Balay if (!BTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 532d5b485abSSatish Balay } 533d5b485abSSatish Balay } 534d5b485abSSatish Balay isz[i] = isz_i; 535d5b485abSSatish Balay } 536*3a40ed3dSBarry Smith PetscFunctionReturn(0); 537d5b485abSSatish Balay } 5385615d1e5SSatish Balay #undef __FUNC__ 539d4bb536fSBarry Smith #define __FUNC__ "MatIncreaseOverlap_MPIBAIJ_Receive" 540d5b485abSSatish Balay /* 541df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 542d5b485abSSatish Balay and return the output 543d5b485abSSatish Balay 544d5b485abSSatish Balay Input: 545d5b485abSSatish Balay C - the matrix 546d5b485abSSatish Balay nrqr - no of messages being processed. 547d5b485abSSatish Balay rbuf - an array of pointers to the recieved requests 548d5b485abSSatish Balay 549d5b485abSSatish Balay Output: 550d5b485abSSatish Balay xdata - array of messages to be sent back 551d5b485abSSatish Balay isz1 - size of each message 552d5b485abSSatish Balay 553d5b485abSSatish Balay For better efficiency perhaps we should malloc seperately each xdata[i], 554d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 555d5b485abSSatish Balay rather then all previous rows as it is now where a single large chunck of 556d5b485abSSatish Balay memory is used. 557d5b485abSSatish Balay 558d5b485abSSatish Balay */ 559df36731bSBarry Smith static int MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,int nrqr,int **rbuf, 560d5b485abSSatish Balay int **xdata, int * isz1) 561d5b485abSSatish Balay { 562df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 563d5b485abSSatish Balay Mat A = c->A, B = c->B; 564df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 565df36731bSBarry Smith int rstart,cstart,*ai, *aj, *bi, *bj, *garray, i, j, k; 566d5b485abSSatish Balay int row,total_sz,ct, ct1, ct2, ct3,mem_estimate, oct2, l, start, end; 567df36731bSBarry Smith int val, max1, max2, rank, Mbs, no_malloc =0, *tmp, new_estimate, ctr; 568bbd702dbSSatish Balay int *rbuf_i,kmax,rbuf_0,ierr; 569bbd702dbSSatish Balay BT xtable; 570d5b485abSSatish Balay 571*3a40ed3dSBarry Smith PetscFunctionBegin; 572d5b485abSSatish Balay rank = c->rank; 573df36731bSBarry Smith Mbs = c->Mbs; 574d5b485abSSatish Balay rstart = c->rstart; 575d5b485abSSatish Balay cstart = c->cstart; 576d5b485abSSatish Balay ai = a->i; 577df36731bSBarry Smith aj = a->j; 578d5b485abSSatish Balay bi = b->i; 579df36731bSBarry Smith bj = b->j; 580d5b485abSSatish Balay garray = c->garray; 581d5b485abSSatish Balay 582d5b485abSSatish Balay 583d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 584d5b485abSSatish Balay rbuf_i = rbuf[i]; 585d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 586d5b485abSSatish Balay ct += rbuf_0; 587d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 588d5b485abSSatish Balay } 589d5b485abSSatish Balay 590df36731bSBarry Smith max1 = ct*(a->nz +b->nz)/c->Mbs; 591d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 592d5b485abSSatish Balay xdata[0] = (int *)PetscMalloc(mem_estimate*sizeof(int)); CHKPTRQ(xdata[0]); 593d5b485abSSatish Balay ++no_malloc; 594bbd702dbSSatish Balay ierr = BTCreate(Mbs,xtable); CHKERRQ(ierr); 595d5b485abSSatish Balay PetscMemzero(isz1,nrqr*sizeof(int)); 596d5b485abSSatish Balay 597d5b485abSSatish Balay ct3 = 0; 598d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 599d5b485abSSatish Balay rbuf_i = rbuf[i]; 600d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 601d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 602d5b485abSSatish Balay ct2 = ct1; 603d5b485abSSatish Balay ct3 += ct1; 604d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 605bbd702dbSSatish Balay BTMemzero(Mbs,xtable); 606d5b485abSSatish Balay oct2 = ct2; 607d5b485abSSatish Balay kmax = rbuf_i[2*j]; 608d5b485abSSatish Balay for (k=0; k<kmax; k++, ct1++) { 609d5b485abSSatish Balay row = rbuf_i[ct1]; 610bbd702dbSSatish Balay if (!BTLookupSet(xtable,row)) { 611d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 612d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 613d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 614d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 615d5b485abSSatish Balay PetscFree(xdata[0]); 616d5b485abSSatish Balay xdata[0] = tmp; 617d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 618d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 619d5b485abSSatish Balay } 620d5b485abSSatish Balay xdata[i][ct2++] = row; 621d5b485abSSatish Balay ct3++; 622d5b485abSSatish Balay } 623d5b485abSSatish Balay } 624d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 625d5b485abSSatish Balay row = xdata[i][k] - rstart; 626d5b485abSSatish Balay start = ai[row]; 627d5b485abSSatish Balay end = ai[row+1]; 628d5b485abSSatish Balay for (l=start; l<end; l++) { 629df36731bSBarry Smith val = aj[l] + cstart; 630bbd702dbSSatish Balay if (!BTLookupSet(xtable,val)) { 631d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 632d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 633d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 634d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 635d5b485abSSatish Balay PetscFree(xdata[0]); 636d5b485abSSatish Balay xdata[0] = tmp; 637d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 638d5b485abSSatish Balay for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 639d5b485abSSatish Balay } 640d5b485abSSatish Balay xdata[i][ct2++] = val; 641d5b485abSSatish Balay ct3++; 642d5b485abSSatish Balay } 643d5b485abSSatish Balay } 644d5b485abSSatish Balay start = bi[row]; 645d5b485abSSatish Balay end = bi[row+1]; 646d5b485abSSatish Balay for (l=start; l<end; l++) { 647df36731bSBarry Smith val = garray[bj[l]]; 648bbd702dbSSatish Balay if (!BTLookupSet(xtable,val)) { 649d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 650d5b485abSSatish Balay new_estimate = (int)(1.5*mem_estimate)+1; 651d5b485abSSatish Balay tmp = (int*) PetscMalloc(new_estimate * sizeof(int));CHKPTRQ(tmp); 652d5b485abSSatish Balay PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int)); 653d5b485abSSatish Balay PetscFree(xdata[0]); 654d5b485abSSatish Balay xdata[0] = tmp; 655d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 656d5b485abSSatish Balay for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 657d5b485abSSatish Balay } 658d5b485abSSatish Balay xdata[i][ct2++] = val; 659d5b485abSSatish Balay ct3++; 660d5b485abSSatish Balay } 661d5b485abSSatish Balay } 662d5b485abSSatish Balay } 663d5b485abSSatish Balay /* Update the header*/ 664d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 665d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 666d5b485abSSatish Balay } 667d5b485abSSatish Balay xdata[i][0] = rbuf_0; 668d5b485abSSatish Balay xdata[i+1] = xdata[i] + ct2; 669d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 670d5b485abSSatish Balay } 671bbd702dbSSatish Balay BTDestroy(xtable); 672df36731bSBarry Smith PLogInfo(0,"MatIncreaseOverlap_MPIBAIJ:[%d] Allocated %d bytes, required %d bytes, no of mallocs = %d\n",rank,mem_estimate, ct3,no_malloc); 673*3a40ed3dSBarry Smith PetscFunctionReturn(0); 674d5b485abSSatish Balay } 675d5b485abSSatish Balay 676cf4f527aSSatish Balay static int MatGetSubMatrices_MPIBAIJ_local(Mat,int,IS *,IS *,MatGetSubMatrixCall,Mat *); 677a2feddc5SSatish Balay 6785615d1e5SSatish Balay #undef __FUNC__ 679d4bb536fSBarry Smith #define __FUNC__ "MatGetSubMatrices_MPIBAIJ" 680df36731bSBarry Smith int MatGetSubMatrices_MPIBAIJ(Mat C,int ismax,IS *isrow,IS *iscol, 681d5b485abSSatish Balay MatGetSubMatrixCall scall,Mat **submat) 682d5b485abSSatish Balay { 68336f4e84dSSatish Balay IS *isrow_new,*iscol_new; 684cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 685cf4f527aSSatish Balay int nmax,nstages_local,nstages,i,pos,max_no,ierr; 686a2feddc5SSatish Balay 687*3a40ed3dSBarry Smith PetscFunctionBegin; 688*3a40ed3dSBarry Smith /* The compression and expansion should be avoided. Does'nt point 689a2feddc5SSatish Balay out errors might change the indices hence buggey */ 690a2feddc5SSatish Balay 69136f4e84dSSatish Balay isrow_new = (IS *)PetscMalloc(2*ismax*sizeof(IS)); CHKPTRQ(isrow_new); 69236f4e84dSSatish Balay iscol_new = isrow_new + ismax; 693e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C, ismax, isrow,isrow_new); CHKERRQ(ierr); 694e757655aSSatish Balay ierr = MatCompressIndicesSorted_MPIBAIJ(C, ismax, iscol,iscol_new); CHKERRQ(ierr); 695cf4f527aSSatish Balay 696cf4f527aSSatish Balay /* Allocate memory to hold all the submatrices */ 697cf4f527aSSatish Balay if (scall != MAT_REUSE_MATRIX) { 698cf4f527aSSatish Balay *submat = (Mat *)PetscMalloc((ismax+1)*sizeof(Mat));CHKPTRQ(*submat); 699cf4f527aSSatish Balay } 700cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 701cf4f527aSSatish Balay nmax = 20*1000000 / (c->Nbs * sizeof(int)); 702cf4f527aSSatish Balay if (nmax == 0) nmax = 1; 703cf4f527aSSatish Balay nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 704cf4f527aSSatish Balay 705cf4f527aSSatish Balay /* Make sure every porcessor loops through the nstages */ 706cf4f527aSSatish Balay MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm); 707cf4f527aSSatish Balay 708cf4f527aSSatish Balay 709cf4f527aSSatish Balay for ( i=0,pos=0; i<nstages; i++ ) { 710cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 711cf4f527aSSatish Balay else if (pos == ismax) max_no = 0; 712cf4f527aSSatish Balay else max_no = ismax-pos; 713cf4f527aSSatish Balay ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_new+pos,iscol_new+pos,scall,*submat+pos); CHKERRQ(ierr); 714cf4f527aSSatish Balay pos += max_no; 715cf4f527aSSatish Balay } 71636f4e84dSSatish Balay 71736f4e84dSSatish Balay for (i=0; i<ismax; i++) { 71836f4e84dSSatish Balay ISDestroy(isrow_new[i]); 71936f4e84dSSatish Balay ISDestroy(iscol_new[i]); 72036f4e84dSSatish Balay } 72136f4e84dSSatish Balay PetscFree(isrow_new); 722*3a40ed3dSBarry Smith PetscFunctionReturn(0); 723a2feddc5SSatish Balay } 724a2feddc5SSatish Balay 725a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 7265615d1e5SSatish Balay #undef __FUNC__ 727d4bb536fSBarry Smith #define __FUNC__ "MatGetSubMatrices_MPIBAIJ_local" 728a2feddc5SSatish Balay static int MatGetSubMatrices_MPIBAIJ_local(Mat C,int ismax,IS *isrow,IS *iscol, 729cf4f527aSSatish Balay MatGetSubMatrixCall scall,Mat *submats) 730a2feddc5SSatish Balay { 731df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *) C->data; 732cf4f527aSSatish Balay Mat A = c->A; 733df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data, *b = (Mat_SeqBAIJ*)c->B->data, *mat; 734d5b485abSSatish Balay int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size; 735a2feddc5SSatish Balay int **sbuf1,**sbuf2, rank, Mbs,i,j,k,l,ct1,ct2,ierr, **rbuf1,row,proc; 73634e6ae18SSatish Balay int nrqs, msz, **ptr,index,*req_size,*ctr,*pa,*tmp,tcol,bsz,nrqr; 737df36731bSBarry Smith int **rbuf3,*req_source,**sbuf_aj, **rbuf2, max1,max2,**rmap; 738d5b485abSSatish Balay int **cmap,**lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 739d5b485abSSatish Balay int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*cmap_i,*lens_i; 74036f4e84dSSatish Balay int *rmap_i,bs=c->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA, *cworkB; 74136f4e84dSSatish Balay int cstart = c->cstart,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 74234e6ae18SSatish Balay int *bmap = c->garray,ctmp,rstart=c->rstart,tag0,tag1,tag2,tag3; 743d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 744d5b485abSSatish Balay MPI_Request *r_waits4,*s_waits3,*s_waits4; 745d5b485abSSatish Balay MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 746d5b485abSSatish Balay MPI_Status *r_status3,*r_status4,*s_status4; 747d5b485abSSatish Balay MPI_Comm comm; 74836f4e84dSSatish Balay Scalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i,*vworkA,*vworkB; 74936f4e84dSSatish Balay Scalar *a_a=a->a,*b_a=b->a; 750d5b485abSSatish Balay 751*3a40ed3dSBarry Smith PetscFunctionBegin; 752d5b485abSSatish Balay comm = C->comm; 75334e6ae18SSatish Balay tag0 = C->tag; 754d5b485abSSatish Balay size = c->size; 755d5b485abSSatish Balay rank = c->rank; 756a2feddc5SSatish Balay Mbs = c->Mbs; 757d5b485abSSatish Balay 75834e6ae18SSatish Balay /* Get some new tags to keep the communication clean */ 75934e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1); CHKERRQ(ierr); 76034e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2); CHKERRQ(ierr); 76134e6ae18SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag3); CHKERRQ(ierr); 76234e6ae18SSatish Balay 763d5b485abSSatish Balay /* Check if the col indices are sorted */ 764d5b485abSSatish Balay for (i=0; i<ismax; i++) { 765d5b485abSSatish Balay ierr = ISSorted(iscol[i],(PetscTruth*)&j); 766e3372554SBarry Smith if (!j) SETERRQ(1,0,"IS is not sorted"); 767d5b485abSSatish Balay } 768d5b485abSSatish Balay 769a2feddc5SSatish Balay len = (2*ismax+1)*(sizeof(int *) + sizeof(int)) + (Mbs+1)*sizeof(int); 770d5b485abSSatish Balay irow = (int **)PetscMalloc(len); CHKPTRQ(irow); 771d5b485abSSatish Balay icol = irow + ismax; 772d5b485abSSatish Balay nrow = (int *) (icol + ismax); 773d5b485abSSatish Balay ncol = nrow + ismax; 774d5b485abSSatish Balay rtable = ncol + ismax; 775d5b485abSSatish Balay 776d5b485abSSatish Balay for (i=0; i<ismax; i++) { 777d5b485abSSatish Balay ierr = ISGetIndices(isrow[i],&irow[i]); CHKERRQ(ierr); 778d5b485abSSatish Balay ierr = ISGetIndices(iscol[i],&icol[i]); CHKERRQ(ierr); 779d5b485abSSatish Balay ierr = ISGetSize(isrow[i],&nrow[i]); CHKERRQ(ierr); 780d5b485abSSatish Balay ierr = ISGetSize(iscol[i],&ncol[i]); CHKERRQ(ierr); 781d5b485abSSatish Balay } 782d5b485abSSatish Balay 783d5b485abSSatish Balay /* Create hash table for the mapping :row -> proc*/ 784d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 785d5b485abSSatish Balay jmax = c->rowners[i+1]; 786d5b485abSSatish Balay for (; j<jmax; j++) { 787d5b485abSSatish Balay rtable[j] = i; 788d5b485abSSatish Balay } 789d5b485abSSatish Balay } 790d5b485abSSatish Balay 791d5b485abSSatish Balay /* evaluate communication - mesg to who, length of mesg, and buffer space 792d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 793d5b485abSSatish Balay w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */ 794d5b485abSSatish Balay w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 795d5b485abSSatish Balay w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 796d5b485abSSatish Balay w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 797d5b485abSSatish Balay PetscMemzero(w1,size*3*sizeof(int)); /* initialise work vector*/ 798d5b485abSSatish Balay for (i=0; i<ismax; i++) { 799d5b485abSSatish Balay PetscMemzero(w4,size*sizeof(int)); /* initialise work vector*/ 800d5b485abSSatish Balay jmax = nrow[i]; 801d5b485abSSatish Balay irow_i = irow[i]; 802d5b485abSSatish Balay for (j=0; j<jmax; j++) { 803d5b485abSSatish Balay row = irow_i[j]; 804d5b485abSSatish Balay proc = rtable[row]; 805d5b485abSSatish Balay w4[proc]++; 806d5b485abSSatish Balay } 807d5b485abSSatish Balay for (j=0; j<size; j++) { 808d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 809d5b485abSSatish Balay } 810d5b485abSSatish Balay } 811d5b485abSSatish Balay 812d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 813e757655aSSatish Balay msz = 0; /* total mesg length for all proc */ 814d5b485abSSatish Balay w1[rank] = 0; /* no mesg sent to intself */ 815d5b485abSSatish Balay w3[rank] = 0; 816d5b485abSSatish Balay for (i=0; i<size; i++) { 817d5b485abSSatish Balay if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 818d5b485abSSatish Balay } 819d5b485abSSatish Balay pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/ 820d5b485abSSatish Balay for (i=0, j=0; i<size; i++) { 821d5b485abSSatish Balay if (w1[i]) { pa[j] = i; j++; } 822d5b485abSSatish Balay } 823d5b485abSSatish Balay 824d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 825d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 826d5b485abSSatish Balay j = pa[i]; 827d5b485abSSatish Balay w1[j] += w2[j] + 2* w3[j]; 828d5b485abSSatish Balay msz += w1[j]; 829d5b485abSSatish Balay } 830d5b485abSSatish Balay /* Do a global reduction to determine how many messages to expect*/ 831d5b485abSSatish Balay { 832d5b485abSSatish Balay int *rw1, *rw2; 833d5b485abSSatish Balay rw1 = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1); 834d5b485abSSatish Balay rw2 = rw1+size; 835d5b485abSSatish Balay MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm); 836d5b485abSSatish Balay bsz = rw1[rank]; 837d5b485abSSatish Balay MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm); 838d5b485abSSatish Balay nrqr = rw2[rank]; 839d5b485abSSatish Balay PetscFree(rw1); 840d5b485abSSatish Balay } 841d5b485abSSatish Balay 842d5b485abSSatish Balay /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 843d5b485abSSatish Balay len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 844d5b485abSSatish Balay rbuf1 = (int**) PetscMalloc(len); CHKPTRQ(rbuf1); 845d5b485abSSatish Balay rbuf1[0] = (int *) (rbuf1 + nrqr); 846d5b485abSSatish Balay for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz; 847d5b485abSSatish Balay 848d5b485abSSatish Balay /* Post the receives */ 8495a197eb2SBarry Smith r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1); 850d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 85134e6ae18SSatish Balay MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i); 852d5b485abSSatish Balay } 853d5b485abSSatish Balay 854d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 855d5b485abSSatish Balay len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 856d5b485abSSatish Balay sbuf1 = (int **)PetscMalloc(len); CHKPTRQ(sbuf1); 857d5b485abSSatish Balay ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 858d5b485abSSatish Balay PetscMemzero(sbuf1,2*size*sizeof(int*)); 859d5b485abSSatish Balay /* allocate memory for outgoing data + buf to receive the first reply */ 860d5b485abSSatish Balay tmp = (int *) (ptr + size); 861d5b485abSSatish Balay ctr = tmp + 2*msz; 862d5b485abSSatish Balay 863d5b485abSSatish Balay { 86436f4e84dSSatish Balay 865d5b485abSSatish Balay int *iptr = tmp,ict = 0; 866d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 867d5b485abSSatish Balay j = pa[i]; 868d5b485abSSatish Balay iptr += ict; 869d5b485abSSatish Balay sbuf1[j] = iptr; 870d5b485abSSatish Balay ict = w1[j]; 871d5b485abSSatish Balay } 872d5b485abSSatish Balay } 873d5b485abSSatish Balay 874d5b485abSSatish Balay /* Form the outgoing messages */ 875d5b485abSSatish Balay /* Initialise the header space */ 876d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 877d5b485abSSatish Balay j = pa[i]; 878d5b485abSSatish Balay sbuf1[j][0] = 0; 879d5b485abSSatish Balay PetscMemzero(sbuf1[j]+1, 2*w3[j]*sizeof(int)); 880d5b485abSSatish Balay ptr[j] = sbuf1[j] + 2*w3[j] + 1; 881d5b485abSSatish Balay } 882d5b485abSSatish Balay 883d5b485abSSatish Balay /* Parse the isrow and copy data into outbuf */ 884d5b485abSSatish Balay for (i=0; i<ismax; i++) { 885d5b485abSSatish Balay PetscMemzero(ctr,size*sizeof(int)); 886d5b485abSSatish Balay irow_i = irow[i]; 887d5b485abSSatish Balay jmax = nrow[i]; 888d5b485abSSatish Balay for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 889d5b485abSSatish Balay row = irow_i[j]; 890d5b485abSSatish Balay proc = rtable[row]; 891d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buf*/ 892d5b485abSSatish Balay ctr[proc]++; 893d5b485abSSatish Balay *ptr[proc] = row; 894d5b485abSSatish Balay ptr[proc]++; 895d5b485abSSatish Balay } 896d5b485abSSatish Balay } 897d5b485abSSatish Balay /* Update the headers for the current IS */ 898d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop too */ 899d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 900d5b485abSSatish Balay sbuf1_j = sbuf1[j]; 901d5b485abSSatish Balay k = ++sbuf1_j[0]; 902d5b485abSSatish Balay sbuf1_j[2*k] = ctr_j; 903d5b485abSSatish Balay sbuf1_j[2*k-1] = i; 904d5b485abSSatish Balay } 905d5b485abSSatish Balay } 906d5b485abSSatish Balay } 907d5b485abSSatish Balay 908d5b485abSSatish Balay /* Now post the sends */ 9095a197eb2SBarry Smith s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1); 910d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 911d5b485abSSatish Balay j = pa[i]; 912d5b485abSSatish Balay /* printf("[%d] Send Req to %d: size %d \n", rank,j, w1[j]); */ 91334e6ae18SSatish Balay MPI_Isend( sbuf1[j], w1[j], MPI_INT, j, tag0, comm, s_waits1+i); 914d5b485abSSatish Balay } 915d5b485abSSatish Balay 916d5b485abSSatish Balay /* Post Recieves to capture the buffer size */ 9175a197eb2SBarry Smith r_waits2 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2); 918d5b485abSSatish Balay rbuf2 = (int**)PetscMalloc((nrqs+1)*sizeof(int *));CHKPTRQ(rbuf2); 919d5b485abSSatish Balay rbuf2[0] = tmp + msz; 920d5b485abSSatish Balay for (i=1; i<nrqs; ++i) { 921d5b485abSSatish Balay j = pa[i]; 922d5b485abSSatish Balay rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 923d5b485abSSatish Balay } 924d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 925d5b485abSSatish Balay j = pa[i]; 92634e6ae18SSatish Balay MPI_Irecv( rbuf2[i], w1[j], MPI_INT, j, tag1, comm, r_waits2+i); 927d5b485abSSatish Balay } 928d5b485abSSatish Balay 929d5b485abSSatish Balay /* Send to other procs the buf size they should allocate */ 930d5b485abSSatish Balay 931d5b485abSSatish Balay 932d5b485abSSatish Balay /* Receive messages*/ 9335a197eb2SBarry Smith s_waits2 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2); 9345a197eb2SBarry Smith r_status1 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1); 935d5b485abSSatish Balay len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*); 936d5b485abSSatish Balay sbuf2 = (int**) PetscMalloc(len);CHKPTRQ(sbuf2); 937d5b485abSSatish Balay req_size = (int *) (sbuf2 + nrqr); 938d5b485abSSatish Balay req_source = req_size + nrqr; 939d5b485abSSatish Balay 940d5b485abSSatish Balay { 941df36731bSBarry Smith Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*) c->A->data, *sB = (Mat_SeqBAIJ*) c->B->data; 94236f4e84dSSatish Balay int *sAi = sA->i, *sBi = sB->i, id; 943d5b485abSSatish Balay int *sbuf2_i; 944d5b485abSSatish Balay 945d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 946d5b485abSSatish Balay MPI_Waitany(nrqr, r_waits1, &index, r_status1+i); 947d5b485abSSatish Balay req_size[index] = 0; 948d5b485abSSatish Balay rbuf1_i = rbuf1[index]; 949d5b485abSSatish Balay start = 2*rbuf1_i[0] + 1; 950d5b485abSSatish Balay MPI_Get_count(r_status1+i,MPI_INT, &end); 951d5b485abSSatish Balay sbuf2[index] = (int *)PetscMalloc(end*sizeof(int));CHKPTRQ(sbuf2[index]); 952d5b485abSSatish Balay sbuf2_i = sbuf2[index]; 953d5b485abSSatish Balay for (j=start; j<end; j++) { 954d5b485abSSatish Balay id = rbuf1_i[j] - rstart; 955d5b485abSSatish Balay ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 956d5b485abSSatish Balay sbuf2_i[j] = ncols; 957d5b485abSSatish Balay req_size[index] += ncols; 958d5b485abSSatish Balay } 959d5b485abSSatish Balay req_source[index] = r_status1[i].MPI_SOURCE; 960d5b485abSSatish Balay /* form the header */ 961d5b485abSSatish Balay sbuf2_i[0] = req_size[index]; 962d5b485abSSatish Balay for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 96334e6ae18SSatish Balay MPI_Isend(sbuf2_i,end,MPI_INT,req_source[index],tag1,comm,s_waits2+i); 964d5b485abSSatish Balay } 965d5b485abSSatish Balay } 966d5b485abSSatish Balay PetscFree(r_status1); PetscFree(r_waits1); 967d5b485abSSatish Balay 968d5b485abSSatish Balay /* recv buffer sizes */ 969d5b485abSSatish Balay /* Receive messages*/ 970d5b485abSSatish Balay 971d5b485abSSatish Balay rbuf3 = (int**)PetscMalloc((nrqs+1)*sizeof(int*)); CHKPTRQ(rbuf3); 972d5b485abSSatish Balay rbuf4 = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar*));CHKPTRQ(rbuf4); 9735a197eb2SBarry Smith r_waits3 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits3); 9745a197eb2SBarry Smith r_waits4 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits4); 9755a197eb2SBarry Smith r_status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2); 976d5b485abSSatish Balay 977d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 978d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits2, &index, r_status2+i); 9795a197eb2SBarry Smith rbuf3[index] = (int *)PetscMalloc(rbuf2[index][0]*sizeof(int));CHKPTRQ(rbuf3[index]); 9805a197eb2SBarry Smith rbuf4[index] = (Scalar *)PetscMalloc(rbuf2[index][0]*bs2*sizeof(Scalar));CHKPTRQ(rbuf4[index]); 981d5b485abSSatish Balay MPI_Irecv(rbuf3[index],rbuf2[index][0], MPI_INT, 98234e6ae18SSatish Balay r_status2[i].MPI_SOURCE, tag2, comm, r_waits3+index); 983a2feddc5SSatish Balay MPI_Irecv(rbuf4[index],rbuf2[index][0]*bs2, MPIU_SCALAR, 98434e6ae18SSatish Balay r_status2[i].MPI_SOURCE, tag3, comm, r_waits4+index); 985d5b485abSSatish Balay } 986d5b485abSSatish Balay PetscFree(r_status2); PetscFree(r_waits2); 987d5b485abSSatish Balay 988d5b485abSSatish Balay /* Wait on sends1 and sends2 */ 9895a197eb2SBarry Smith s_status1 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1); 9905a197eb2SBarry Smith s_status2 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status2); 991d5b485abSSatish Balay 992d5b485abSSatish Balay MPI_Waitall(nrqs,s_waits1,s_status1); 993d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits2,s_status2); 994d5b485abSSatish Balay PetscFree(s_status1); PetscFree(s_status2); 995d5b485abSSatish Balay PetscFree(s_waits1); PetscFree(s_waits2); 996d5b485abSSatish Balay 997d5b485abSSatish Balay /* Now allocate buffers for a->j, and send them off */ 998d5b485abSSatish Balay sbuf_aj = (int **)PetscMalloc((nrqr+1)*sizeof(int *));CHKPTRQ(sbuf_aj); 999d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1000d5b485abSSatish Balay sbuf_aj[0] = (int*) PetscMalloc((j+1)*sizeof(int)); CHKPTRQ(sbuf_aj[0]); 1001d5b485abSSatish Balay for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1002d5b485abSSatish Balay 10035a197eb2SBarry Smith s_waits3 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits3); 1004d5b485abSSatish Balay { 1005d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1006d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1007d5b485abSSatish Balay sbuf_aj_i = sbuf_aj[i]; 1008d5b485abSSatish Balay ct1 = 2*rbuf1_i[0] + 1; 1009d5b485abSSatish Balay ct2 = 0; 1010d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1011d5b485abSSatish Balay kmax = rbuf1[i][2*j]; 1012d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1013e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1014d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1015d5b485abSSatish Balay ncols = nzA + nzB; 1016d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1017d5b485abSSatish Balay 1018d5b485abSSatish Balay /* load the column indices for this row into cols*/ 1019d5b485abSSatish Balay cols = sbuf_aj_i + ct2; 1020d5b485abSSatish Balay for (l=0; l<nzB; l++) { 1021d5b485abSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1022d5b485abSSatish Balay else break; 1023d5b485abSSatish Balay } 1024d5b485abSSatish Balay imark = l; 1025d5b485abSSatish Balay for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1026d5b485abSSatish Balay for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1027d5b485abSSatish Balay ct2 += ncols; 1028d5b485abSSatish Balay } 1029d5b485abSSatish Balay } 103034e6ae18SSatish Balay MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i); 1031d5b485abSSatish Balay } 1032d5b485abSSatish Balay } 10335a197eb2SBarry Smith r_status3 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status3); 10345a197eb2SBarry Smith s_status3 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status3); 1035d5b485abSSatish Balay 1036d5b485abSSatish Balay /* Allocate buffers for a->a, and send them off */ 1037d5b485abSSatish Balay sbuf_aa = (Scalar **)PetscMalloc((nrqr+1)*sizeof(Scalar *));CHKPTRQ(sbuf_aa); 1038d5b485abSSatish Balay for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1039a2feddc5SSatish Balay sbuf_aa[0] = (Scalar*) PetscMalloc((j+1)*bs2*sizeof(Scalar));CHKPTRQ(sbuf_aa[0]); 1040a2feddc5SSatish Balay for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1041d5b485abSSatish Balay 10425a197eb2SBarry Smith s_waits4 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits4); 1043d5b485abSSatish Balay { 1044d5b485abSSatish Balay for (i=0; i<nrqr; i++) { 1045d5b485abSSatish Balay rbuf1_i = rbuf1[i]; 1046d5b485abSSatish Balay sbuf_aa_i = sbuf_aa[i]; 1047d5b485abSSatish Balay ct1 = 2*rbuf1_i[0]+1; 1048d5b485abSSatish Balay ct2 = 0; 1049d5b485abSSatish Balay for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1050d5b485abSSatish Balay kmax = rbuf1_i[2*j]; 1051d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 1052e8e0db45SSatish Balay row = rbuf1_i[ct1] - rstart; 1053d5b485abSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1054d5b485abSSatish Balay ncols = nzA + nzB; 1055d5b485abSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1056a2feddc5SSatish Balay vworkA = a_a + a_i[row]*bs2; vworkB = b_a + b_i[row]*bs2; 1057d5b485abSSatish Balay 1058d5b485abSSatish Balay /* load the column values for this row into vals*/ 10595b83ace0SSatish Balay vals = sbuf_aa_i+ct2*bs2; 1060d5b485abSSatish Balay for (l=0; l<nzB; l++) { 1061*3a40ed3dSBarry Smith if ((bmap[cworkB[l]]) < cstart) { 1062a2feddc5SSatish Balay PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(Scalar)); 1063*3a40ed3dSBarry Smith } 1064d5b485abSSatish Balay else break; 1065d5b485abSSatish Balay } 1066d5b485abSSatish Balay imark = l; 1067*3a40ed3dSBarry Smith for (l=0; l<nzA; l++) { 1068a2feddc5SSatish Balay PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(Scalar)); 1069*3a40ed3dSBarry Smith } 1070*3a40ed3dSBarry Smith for (l=imark; l<nzB; l++) { 1071a2feddc5SSatish Balay PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(Scalar)); 1072*3a40ed3dSBarry Smith } 1073d5b485abSSatish Balay ct2 += ncols; 1074d5b485abSSatish Balay } 1075d5b485abSSatish Balay } 107634e6ae18SSatish Balay MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i); 1077d5b485abSSatish Balay } 1078d5b485abSSatish Balay } 10795a197eb2SBarry Smith r_status4 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status4); 10805a197eb2SBarry Smith s_status4 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status4); 1081d5b485abSSatish Balay PetscFree(rbuf1); 1082d5b485abSSatish Balay 1083d5b485abSSatish Balay /* Form the matrix */ 1084d5b485abSSatish Balay /* create col map */ 1085d5b485abSSatish Balay { 1086d5b485abSSatish Balay int *icol_i; 1087d5b485abSSatish Balay 1088a2feddc5SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->Nbs*sizeof(int); 1089d5b485abSSatish Balay cmap = (int **)PetscMalloc(len); CHKPTRQ(cmap); 1090d5b485abSSatish Balay cmap[0] = (int *)(cmap + ismax); 1091a2feddc5SSatish Balay PetscMemzero(cmap[0],(1+ismax*c->Nbs)*sizeof(int)); 1092a2feddc5SSatish Balay for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + c->Nbs; } 1093d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1094d5b485abSSatish Balay jmax = ncol[i]; 1095d5b485abSSatish Balay icol_i = icol[i]; 1096d5b485abSSatish Balay cmap_i = cmap[i]; 1097d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1098d5b485abSSatish Balay cmap_i[icol_i[j]] = j+1; 1099d5b485abSSatish Balay } 1100d5b485abSSatish Balay } 1101d5b485abSSatish Balay } 1102d5b485abSSatish Balay 1103d5b485abSSatish Balay 1104d5b485abSSatish Balay /* Create lens which is required for MatCreate... */ 1105d5b485abSSatish Balay for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1106d5b485abSSatish Balay len = (1+ismax)*sizeof(int *) + j*sizeof(int); 1107d5b485abSSatish Balay lens = (int **)PetscMalloc(len); CHKPTRQ(lens); 1108d5b485abSSatish Balay lens[0] = (int *)(lens + ismax); 1109d5b485abSSatish Balay PetscMemzero(lens[0], j*sizeof(int)); 1110d5b485abSSatish Balay for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1111d5b485abSSatish Balay 1112d5b485abSSatish Balay /* Update lens from local data */ 1113d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1114d5b485abSSatish Balay jmax = nrow[i]; 1115d5b485abSSatish Balay cmap_i = cmap[i]; 1116d5b485abSSatish Balay irow_i = irow[i]; 1117d5b485abSSatish Balay lens_i = lens[i]; 1118d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1119d5b485abSSatish Balay row = irow_i[j]; 1120d5b485abSSatish Balay proc = rtable[row]; 1121d5b485abSSatish Balay if (proc == rank) { 112236f4e84dSSatish Balay /* Get indices from matA and then from matB */ 112336f4e84dSSatish Balay row = row - rstart; 112436f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 112536f4e84dSSatish Balay cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 112636f4e84dSSatish Balay for (k=0; k<nzA; k++) 112736f4e84dSSatish Balay if (cmap_i[cstart + cworkA[k]]) { lens_i[j]++;} 112836f4e84dSSatish Balay for (k=0; k<nzB; k++) 112936f4e84dSSatish Balay if (cmap_i[bmap[cworkB[k]]]) { lens_i[j]++;} 1130d5b485abSSatish Balay } 1131d5b485abSSatish Balay } 1132a2feddc5SSatish Balay } 1133d5b485abSSatish Balay 1134d5b485abSSatish Balay /* Create row map*/ 1135a2feddc5SSatish Balay len = (1+ismax)*sizeof(int *) + ismax*c->Mbs*sizeof(int); 1136d5b485abSSatish Balay rmap = (int **)PetscMalloc(len); CHKPTRQ(rmap); 1137d5b485abSSatish Balay rmap[0] = (int *)(rmap + ismax); 1138a2feddc5SSatish Balay PetscMemzero(rmap[0],ismax*c->Mbs*sizeof(int)); 1139a2feddc5SSatish Balay for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + c->Mbs;} 1140d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1141d5b485abSSatish Balay rmap_i = rmap[i]; 1142d5b485abSSatish Balay irow_i = irow[i]; 1143d5b485abSSatish Balay jmax = nrow[i]; 1144d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1145d5b485abSSatish Balay rmap_i[irow_i[j]] = j; 1146d5b485abSSatish Balay } 1147d5b485abSSatish Balay } 1148d5b485abSSatish Balay 1149d5b485abSSatish Balay /* Update lens from offproc data */ 1150d5b485abSSatish Balay { 1151d5b485abSSatish Balay int *rbuf2_i, *rbuf3_i, *sbuf1_i; 1152d5b485abSSatish Balay 1153d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1154d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits3, &i, r_status3+tmp2); 1155d5b485abSSatish Balay index = pa[i]; 1156d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1157d5b485abSSatish Balay jmax = sbuf1_i[0]; 1158d5b485abSSatish Balay ct1 = 2*jmax+1; 1159d5b485abSSatish Balay ct2 = 0; 1160d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1161d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1162d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1163d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1164d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1165d5b485abSSatish Balay lens_i = lens[is_no]; 1166d5b485abSSatish Balay cmap_i = cmap[is_no]; 1167d5b485abSSatish Balay rmap_i = rmap[is_no]; 1168d5b485abSSatish Balay for (k=0; k<max1; k++,ct1++) { 1169d5b485abSSatish Balay row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1170d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1171d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1172d5b485abSSatish Balay if (cmap_i[rbuf3_i[ct2]]) { 1173d5b485abSSatish Balay lens_i[row]++; 1174d5b485abSSatish Balay } 1175d5b485abSSatish Balay } 1176d5b485abSSatish Balay } 1177d5b485abSSatish Balay } 1178d5b485abSSatish Balay } 1179d5b485abSSatish Balay } 1180d5b485abSSatish Balay PetscFree(r_status3); PetscFree(r_waits3); 1181d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits3,s_status3); 1182d5b485abSSatish Balay PetscFree(s_status3); PetscFree(s_waits3); 1183d5b485abSSatish Balay 1184d5b485abSSatish Balay /* Create the submatrices */ 1185d5b485abSSatish Balay if (scall == MAT_REUSE_MATRIX) { 1186d5b485abSSatish Balay /* 1187d5b485abSSatish Balay Assumes new rows are same length as the old rows, hence bug! 1188d5b485abSSatish Balay */ 1189d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1190df36731bSBarry Smith mat = (Mat_SeqBAIJ *)(submats[i]->data); 119136f4e84dSSatish Balay if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || mat->bs != bs)) { 1192e3372554SBarry Smith SETERRQ(1,0,"Cannot reuse matrix. wrong size"); 1193d5b485abSSatish Balay } 119436f4e84dSSatish Balay if (PetscMemcmp(mat->ilen,lens[i], mat->mbs *sizeof(int))) { 1195e3372554SBarry Smith SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros"); 1196d5b485abSSatish Balay } 1197d5b485abSSatish Balay /* Initial matrix as if empty */ 119836f4e84dSSatish Balay PetscMemzero(mat->ilen,mat->mbs*sizeof(int)); 1199ce60de66SLois Curfman McInnes submats[i]->factor = C->factor; 1200d5b485abSSatish Balay } 1201d5b485abSSatish Balay } 1202d5b485abSSatish Balay else { 1203cf4f527aSSatish Balay /* *submat = submats = (Mat *)PetscMalloc(ismax*sizeof(Mat)); CHKPTRQ(submats); */ 1204d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1205029af93fSBarry Smith ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,a->bs,nrow[i]*bs,ncol[i]*bs,0,lens[i],submats+i); 12063b2fbd54SBarry Smith CHKERRQ(ierr); 1207d5b485abSSatish Balay } 1208d5b485abSSatish Balay } 1209d5b485abSSatish Balay 1210d5b485abSSatish Balay /* Assemble the matrices */ 1211d5b485abSSatish Balay /* First assemble the local rows */ 1212d5b485abSSatish Balay { 121336f4e84dSSatish Balay int ilen_row,*imat_ilen, *imat_j, *imat_i; 1214d5b485abSSatish Balay Scalar *imat_a; 1215d5b485abSSatish Balay 1216d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1217df36731bSBarry Smith mat = (Mat_SeqBAIJ *) submats[i]->data; 1218d5b485abSSatish Balay imat_ilen = mat->ilen; 1219d5b485abSSatish Balay imat_j = mat->j; 1220d5b485abSSatish Balay imat_i = mat->i; 1221d5b485abSSatish Balay imat_a = mat->a; 1222d5b485abSSatish Balay cmap_i = cmap[i]; 1223d5b485abSSatish Balay rmap_i = rmap[i]; 1224d5b485abSSatish Balay irow_i = irow[i]; 1225d5b485abSSatish Balay jmax = nrow[i]; 1226d5b485abSSatish Balay for (j=0; j<jmax; j++) { 1227d5b485abSSatish Balay row = irow_i[j]; 1228d5b485abSSatish Balay proc = rtable[row]; 1229d5b485abSSatish Balay if (proc == rank) { 123036f4e84dSSatish Balay row = row - rstart; 123136f4e84dSSatish Balay nzA = a_i[row+1] - a_i[row]; 123236f4e84dSSatish Balay nzB = b_i[row+1] - b_i[row]; 123336f4e84dSSatish Balay cworkA = a_j + a_i[row]; 123436f4e84dSSatish Balay cworkB = b_j + b_i[row]; 123536f4e84dSSatish Balay vworkA = a_a + a_i[row]*bs2; 123636f4e84dSSatish Balay vworkB = b_a + b_i[row]*bs2; 123736f4e84dSSatish Balay 123836f4e84dSSatish Balay row = rmap_i[row + rstart]; 1239df36731bSBarry Smith mat_i = imat_i[row]; 124036f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1241d5b485abSSatish Balay mat_j = imat_j + mat_i; 124236f4e84dSSatish Balay ilen_row = imat_ilen[row]; 124336f4e84dSSatish Balay 124436f4e84dSSatish Balay /* load the column indices for this row into cols*/ 124536f4e84dSSatish Balay for (l=0; l<nzB; l++) { 124636f4e84dSSatish Balay if ((ctmp = bmap[cworkB[l]]) < cstart) { 124736f4e84dSSatish Balay if ((tcol = cmap_i[ctmp])) { 1248df36731bSBarry Smith *mat_j++ = tcol - 1; 124936f4e84dSSatish Balay PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 1250d5b485abSSatish Balay ilen_row++; 1251d5b485abSSatish Balay } 1252d5b485abSSatish Balay } 125336f4e84dSSatish Balay else break; 125436f4e84dSSatish Balay } 125536f4e84dSSatish Balay imark = l; 125636f4e84dSSatish Balay for (l=0; l<nzA; l++) { 125736f4e84dSSatish Balay if ((tcol = cmap_i[cstart + cworkA[l]])) { 125836f4e84dSSatish Balay *mat_j++ = tcol - 1; 125936f4e84dSSatish Balay PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 126036f4e84dSSatish Balay ilen_row++; 126136f4e84dSSatish Balay } 126236f4e84dSSatish Balay } 126336f4e84dSSatish Balay for (l=imark; l<nzB; l++) { 126436f4e84dSSatish Balay if ((tcol = cmap_i[bmap[cworkB[l]]])) { 126536f4e84dSSatish Balay *mat_j++ = tcol - 1; 126636f4e84dSSatish Balay PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 126736f4e84dSSatish Balay ilen_row++; 126836f4e84dSSatish Balay } 126936f4e84dSSatish Balay } 1270d5b485abSSatish Balay imat_ilen[row] = ilen_row; 1271d5b485abSSatish Balay } 1272d5b485abSSatish Balay } 127336f4e84dSSatish Balay 1274d5b485abSSatish Balay } 1275d5b485abSSatish Balay } 1276d5b485abSSatish Balay 1277d5b485abSSatish Balay /* Now assemble the off proc rows*/ 1278d5b485abSSatish Balay { 1279d5b485abSSatish Balay int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1280d5b485abSSatish Balay int *imat_j,*imat_i; 1281d5b485abSSatish Balay Scalar *imat_a,*rbuf4_i; 1282d5b485abSSatish Balay 1283d5b485abSSatish Balay for (tmp2=0; tmp2<nrqs; tmp2++) { 1284d5b485abSSatish Balay MPI_Waitany(nrqs, r_waits4, &i, r_status4+tmp2); 1285d5b485abSSatish Balay index = pa[i]; 1286d5b485abSSatish Balay sbuf1_i = sbuf1[index]; 1287d5b485abSSatish Balay jmax = sbuf1_i[0]; 1288d5b485abSSatish Balay ct1 = 2*jmax + 1; 1289d5b485abSSatish Balay ct2 = 0; 1290d5b485abSSatish Balay rbuf2_i = rbuf2[i]; 1291d5b485abSSatish Balay rbuf3_i = rbuf3[i]; 1292d5b485abSSatish Balay rbuf4_i = rbuf4[i]; 1293d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 1294d5b485abSSatish Balay is_no = sbuf1_i[2*j-1]; 1295d5b485abSSatish Balay rmap_i = rmap[is_no]; 1296d5b485abSSatish Balay cmap_i = cmap[is_no]; 1297df36731bSBarry Smith mat = (Mat_SeqBAIJ *) submats[is_no]->data; 1298d5b485abSSatish Balay imat_ilen = mat->ilen; 1299d5b485abSSatish Balay imat_j = mat->j; 1300d5b485abSSatish Balay imat_i = mat->i; 1301d5b485abSSatish Balay imat_a = mat->a; 1302d5b485abSSatish Balay max1 = sbuf1_i[2*j]; 1303d5b485abSSatish Balay for (k=0; k<max1; k++, ct1++) { 1304d5b485abSSatish Balay row = sbuf1_i[ct1]; 1305d5b485abSSatish Balay row = rmap_i[row]; 1306d5b485abSSatish Balay ilen = imat_ilen[row]; 1307df36731bSBarry Smith mat_i = imat_i[row]; 130836f4e84dSSatish Balay mat_a = imat_a + mat_i*bs2; 1309d5b485abSSatish Balay mat_j = imat_j + mat_i; 1310d5b485abSSatish Balay max2 = rbuf2_i[ct1]; 1311d5b485abSSatish Balay for (l=0; l<max2; l++,ct2++) { 1312d5b485abSSatish Balay if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1313df36731bSBarry Smith *mat_j++ = tcol - 1; 131436f4e84dSSatish Balay /* *mat_a++ = rbuf4_i[ct2]; */ 131536f4e84dSSatish Balay PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(Scalar)); mat_a += bs2; 1316d5b485abSSatish Balay ilen++; 1317d5b485abSSatish Balay } 1318d5b485abSSatish Balay } 1319d5b485abSSatish Balay imat_ilen[row] = ilen; 1320d5b485abSSatish Balay } 1321d5b485abSSatish Balay } 1322d5b485abSSatish Balay } 1323d5b485abSSatish Balay } 1324d5b485abSSatish Balay PetscFree(r_status4); PetscFree(r_waits4); 1325d5b485abSSatish Balay MPI_Waitall(nrqr,s_waits4,s_status4); 1326d5b485abSSatish Balay PetscFree(s_waits4); PetscFree(s_status4); 1327d5b485abSSatish Balay 1328d5b485abSSatish Balay /* Restore the indices */ 1329d5b485abSSatish Balay for (i=0; i<ismax; i++) { 1330d5b485abSSatish Balay ierr = ISRestoreIndices(isrow[i], irow+i); CHKERRQ(ierr); 1331d5b485abSSatish Balay ierr = ISRestoreIndices(iscol[i], icol+i); CHKERRQ(ierr); 1332d5b485abSSatish Balay } 1333d5b485abSSatish Balay 1334d5b485abSSatish Balay /* Destroy allocated memory */ 1335d5b485abSSatish Balay PetscFree(irow); 1336d5b485abSSatish Balay PetscFree(w1); 1337d5b485abSSatish Balay PetscFree(pa); 1338d5b485abSSatish Balay 1339d5b485abSSatish Balay PetscFree(sbuf1); 1340d5b485abSSatish Balay PetscFree(rbuf2); 1341d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 1342d5b485abSSatish Balay PetscFree(sbuf2[i]); 1343d5b485abSSatish Balay } 1344d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 1345d5b485abSSatish Balay PetscFree(rbuf3[i]); 1346d5b485abSSatish Balay PetscFree(rbuf4[i]); 1347d5b485abSSatish Balay } 1348d5b485abSSatish Balay 1349d5b485abSSatish Balay PetscFree(sbuf2); 1350d5b485abSSatish Balay PetscFree(rbuf3); 1351d5b485abSSatish Balay PetscFree(rbuf4 ); 1352d5b485abSSatish Balay PetscFree(sbuf_aj[0]); 1353d5b485abSSatish Balay PetscFree(sbuf_aj); 1354d5b485abSSatish Balay PetscFree(sbuf_aa[0]); 1355d5b485abSSatish Balay PetscFree(sbuf_aa); 1356d5b485abSSatish Balay 1357d5b485abSSatish Balay PetscFree(cmap); 1358d5b485abSSatish Balay PetscFree(rmap); 1359d5b485abSSatish Balay PetscFree(lens); 1360d5b485abSSatish Balay 1361d5b485abSSatish Balay for (i=0; i<ismax; i++) { 136236f4e84dSSatish Balay ierr = MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 136336f4e84dSSatish Balay ierr = MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1364d5b485abSSatish Balay } 136534e6ae18SSatish Balay 136634e6ae18SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag3); CHKERRQ(ierr); 136734e6ae18SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag2); CHKERRQ(ierr); 136834e6ae18SSatish Balay ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1); CHKERRQ(ierr); 136934e6ae18SSatish Balay 1370*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1371d5b485abSSatish Balay } 1372