123bdbc58SBarry Smith 2d5b485abSSatish Balay /* 3d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 4d5b485abSSatish Balay and to find submatrices that were shared across processors. 5d5b485abSSatish Balay */ 6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 7c6db04a5SJed Brown #include <petscbt.h> 8d5b485abSSatish Balay 9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char**,PetscInt*,PetscInt**); 10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt**,PetscInt**,PetscInt*); 1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 13d5b485abSSatish Balay 14b24ad042SBarry Smith PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 15d5b485abSSatish Balay { 166849ba73SBarry Smith PetscErrorCode ierr; 17d0f46423SBarry Smith PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 1836f4e84dSSatish Balay IS *is_new; 1936f4e84dSSatish Balay 203a40ed3dSBarry Smith PetscFunctionBegin; 21785e854fSJed Brown ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr); 22df36731bSBarry Smith /* Convert the indices into block format */ 2305d8c843SHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 24546078acSJacob Faibussowitsch if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 25d5b485abSSatish Balay for (i=0; i<ov; ++i) { 2636f4e84dSSatish Balay ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 27d5b485abSSatish Balay } 286bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 2905d8c843SHong Zhang ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr); 306bf464f9SBarry Smith for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 31606d414cSSatish Balay ierr = PetscFree(is_new);CHKERRQ(ierr); 323a40ed3dSBarry Smith PetscFunctionReturn(0); 33d5b485abSSatish Balay } 34d5b485abSSatish Balay 35d5b485abSSatish Balay /* 36d5b485abSSatish Balay Sample message format: 37d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 380e9b0e7eSHong Zhang to index sets is[1], is[5] 39d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 40d5b485abSSatish Balay ----------- 41d5b485abSSatish Balay mesg [1] = 1 => is[1] 42d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 43d5b485abSSatish Balay ----------- 44d5b485abSSatish Balay mesg [5] = 5 => is[5] 45d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 46d5b485abSSatish Balay ----------- 47d5b485abSSatish Balay mesg [7] 480e9b0e7eSHong Zhang mesg [n] data(is[1]) 49d5b485abSSatish Balay ----------- 50d5b485abSSatish Balay mesg[n+1] 51d5b485abSSatish Balay mesg[m] data(is[5]) 52d5b485abSSatish Balay ----------- 53d5b485abSSatish Balay 54d5b485abSSatish Balay Notes: 55d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 562d4ee042Sprj- nrqr - no of requests received (which have to be or which have been processed) 57d5b485abSSatish Balay */ 58db41cccfSHong Zhang PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 59d5b485abSSatish Balay { 60df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 615d0c19d7SBarry Smith const PetscInt **idx,*idx_i; 6224c049a4SHong Zhang PetscInt *n,*w3,*w4,**data,len; 636849ba73SBarry Smith PetscErrorCode ierr; 64b24ad042SBarry Smith PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 65131c27b5Sprj- PetscInt Mbs,i,j,k,**rbuf,row,nrqs,msz,**outdat,**ptr; 668f9f447aSHong Zhang PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p; 67131c27b5Sprj- PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2,proc=-1; 68f1af5d2fSBarry Smith PetscBT *table; 69749130a8SStefano Zampini MPI_Comm comm,*iscomms; 70d5b485abSSatish Balay MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 718f9f447aSHong Zhang char *t_p; 72d5b485abSSatish Balay 733a40ed3dSBarry Smith PetscFunctionBegin; 74ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 75d5b485abSSatish Balay size = c->size; 76d5b485abSSatish Balay rank = c->rank; 77df36731bSBarry Smith Mbs = c->Mbs; 78d5b485abSSatish Balay 79c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 80c7dd2924SSatish Balay ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 81c7dd2924SSatish Balay 82f489ac74SBarry Smith ierr = PetscMalloc2(imax+1,(PetscInt***)&idx,imax,&n);CHKERRQ(ierr); 8324c049a4SHong Zhang 84d5b485abSSatish Balay for (i=0; i<imax; i++) { 85d5b485abSSatish Balay ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 86b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 87d5b485abSSatish Balay } 88d5b485abSSatish Balay 89d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 90d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 91884e879aSBarry Smith ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 92d5b485abSSatish Balay for (i=0; i<imax; i++) { 93580bdb30SBarry Smith ierr = PetscArrayzero(w4,size);CHKERRQ(ierr); /* initialise work vector*/ 94d5b485abSSatish Balay idx_i = idx[i]; 95d5b485abSSatish Balay len = n[i]; 96d5b485abSSatish Balay for (j=0; j<len; j++) { 97d5b485abSSatish Balay row = idx_i[j]; 98f23aa3ddSBarry Smith if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 99ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 100d5b485abSSatish Balay w4[proc]++; 101d5b485abSSatish Balay } 102d5b485abSSatish Balay for (j=0; j<size; j++) { 103d5b485abSSatish Balay if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 104d5b485abSSatish Balay } 105d5b485abSSatish Balay } 106d5b485abSSatish Balay 107d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 108d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1090e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 110d5b485abSSatish Balay w3[rank] = 0; 111d5b485abSSatish Balay for (i=0; i<size; i++) { 112d5b485abSSatish Balay if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 113d5b485abSSatish Balay } 114d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 115175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&pa);CHKERRQ(ierr); 116d5b485abSSatish Balay for (i=0,j=0; i<size; i++) { 117d5b485abSSatish Balay if (w1[i]) {pa[j] = i; j++;} 118d5b485abSSatish Balay } 119d5b485abSSatish Balay 120d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 121d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 122d5b485abSSatish Balay j = pa[i]; 123d5b485abSSatish Balay w1[j] += w2[j] + 2*w3[j]; 124d5b485abSSatish Balay msz += w1[j]; 125d5b485abSSatish Balay } 126d5b485abSSatish Balay 127c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 128a96d083dSSatish Balay ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 129c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 130d5b485abSSatish Balay 131c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 132c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 133d5b485abSSatish Balay 134d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 135dcca6d9dSJed Brown ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr); 136580bdb30SBarry Smith ierr = PetscArrayzero(outdat,size);CHKERRQ(ierr); 137580bdb30SBarry Smith ierr = PetscArrayzero(ptr,size);CHKERRQ(ierr); 138d5b485abSSatish Balay { 139b24ad042SBarry Smith PetscInt *iptr = tmp,ict = 0; 140d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 141d5b485abSSatish Balay j = pa[i]; 142d5b485abSSatish Balay iptr += ict; 143d5b485abSSatish Balay outdat[j] = iptr; 144d5b485abSSatish Balay ict = w1[j]; 145d5b485abSSatish Balay } 146d5b485abSSatish Balay } 147d5b485abSSatish Balay 148d5b485abSSatish Balay /* Form the outgoing messages */ 149d5b485abSSatish Balay /*plug in the headers*/ 150d5b485abSSatish Balay for (i=0; i<nrqs; i++) { 151d5b485abSSatish Balay j = pa[i]; 152d5b485abSSatish Balay outdat[j][0] = 0; 153580bdb30SBarry Smith ierr = PetscArrayzero(outdat[j]+1,2*w3[j]);CHKERRQ(ierr); 154d5b485abSSatish Balay ptr[j] = outdat[j] + 2*w3[j] + 1; 155d5b485abSSatish Balay } 156d5b485abSSatish Balay 157d5b485abSSatish Balay /* Memory for doing local proc's work*/ 158d5b485abSSatish Balay { 1591795a4d1SJed Brown ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr); 160d5b485abSSatish Balay 161bbd702dbSSatish Balay for (i=0; i<imax; i++) { 162b6410449SSatish Balay table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 163bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 164d5b485abSSatish Balay } 165d5b485abSSatish Balay } 166d5b485abSSatish Balay 167d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 168d5b485abSSatish Balay { 169b24ad042SBarry Smith PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 170f1af5d2fSBarry Smith PetscBT table_i; 171d5b485abSSatish Balay 172d5b485abSSatish Balay for (i=0; i<imax; i++) { 173580bdb30SBarry Smith ierr = PetscArrayzero(ctr,size);CHKERRQ(ierr); 174d5b485abSSatish Balay n_i = n[i]; 175d5b485abSSatish Balay table_i = table[i]; 176d5b485abSSatish Balay idx_i = idx[i]; 177d5b485abSSatish Balay data_i = data[i]; 178d5b485abSSatish Balay isz_i = isz[i]; 179d5b485abSSatish Balay for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 180d5b485abSSatish Balay row = idx_i[j]; 181ca31476aSJed Brown ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 182d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 183d5b485abSSatish Balay ctr[proc]++; 184d5b485abSSatish Balay *ptr[proc] = row; 185d5b485abSSatish Balay ptr[proc]++; 186d6b45a43SBarry Smith } else { /* Update the local table */ 18726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 188d5b485abSSatish Balay } 189d5b485abSSatish Balay } 190d5b485abSSatish Balay /* Update the headers for the current IS */ 191d5b485abSSatish Balay for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 192d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 193d5b485abSSatish Balay outdat_j = outdat[j]; 194d5b485abSSatish Balay k = ++outdat_j[0]; 195d5b485abSSatish Balay outdat_j[2*k] = ctr_j; 196d5b485abSSatish Balay outdat_j[2*k-1] = i; 197d5b485abSSatish Balay } 198d5b485abSSatish Balay } 199d5b485abSSatish Balay isz[i] = isz_i; 200d5b485abSSatish Balay } 201d5b485abSSatish Balay } 202d5b485abSSatish Balay 203d5b485abSSatish Balay /* Now post the sends */ 204175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&s_waits1);CHKERRQ(ierr); 205d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 206d5b485abSSatish Balay j = pa[i]; 207ffc4695bSBarry Smith ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRMPI(ierr); 208d5b485abSSatish Balay } 209d5b485abSSatish Balay 210d5b485abSSatish Balay /* No longer need the original indices*/ 211d5b485abSSatish Balay for (i=0; i<imax; ++i) { 212d5b485abSSatish Balay ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 213d5b485abSSatish Balay } 2142cd485c2SBarry Smith ierr = PetscFree2(*(PetscInt***)&idx,n);CHKERRQ(ierr); 215d5b485abSSatish Balay 216749130a8SStefano Zampini ierr = PetscMalloc1(imax,&iscomms);CHKERRQ(ierr); 217d5b485abSSatish Balay for (i=0; i<imax; ++i) { 218749130a8SStefano Zampini ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomms[i],NULL);CHKERRQ(ierr); 2196bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 220d5b485abSSatish Balay } 221d5b485abSSatish Balay 222d5b485abSSatish Balay /* Do Local work*/ 223df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 224d5b485abSSatish Balay 225d5b485abSSatish Balay /* Receive messages*/ 226175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqr,r_waits1,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 227175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqs,s_waits1,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 228d5b485abSSatish Balay 229d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 23023bdbc58SBarry Smith ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 23123bdbc58SBarry Smith ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 232d5b485abSSatish Balay 233175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&xdata);CHKERRQ(ierr); 234175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&isz1);CHKERRQ(ierr); 235df36731bSBarry Smith ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 236175c2bc5SJunchao Zhang if (rbuf) { 237c05d87d6SBarry Smith ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 238606d414cSSatish Balay ierr = PetscFree(rbuf);CHKERRQ(ierr); 239175c2bc5SJunchao Zhang } 240d5b485abSSatish Balay 241d5b485abSSatish Balay /* Send the data back*/ 242d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 243d5b485abSSatish Balay { 244b24ad042SBarry Smith PetscMPIInt *rw1; 245d5b485abSSatish Balay 246884e879aSBarry Smith ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr); 247c7dd2924SSatish Balay 248d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 249175c2bc5SJunchao Zhang proc = onodes1[i]; 250d5b485abSSatish Balay rw1[proc] = isz1[i]; 251d5b485abSSatish Balay } 252d5b485abSSatish Balay 253c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 254c7dd2924SSatish Balay ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 25503512dc5SSatish Balay ierr = PetscFree(rw1);CHKERRQ(ierr); 256c7dd2924SSatish Balay } 257c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 258c7dd2924SSatish Balay ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 259d5b485abSSatish Balay 260d5b485abSSatish Balay /* Now post the sends */ 261175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&s_waits2);CHKERRQ(ierr); 262d5b485abSSatish Balay for (i=0; i<nrqr; ++i) { 263175c2bc5SJunchao Zhang j = onodes1[i]; 264ffc4695bSBarry Smith ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRMPI(ierr); 265d5b485abSSatish Balay } 266d5b485abSSatish Balay 267175c2bc5SJunchao Zhang ierr = PetscFree(onodes1);CHKERRQ(ierr); 268175c2bc5SJunchao Zhang ierr = PetscFree(olengths1);CHKERRQ(ierr); 269175c2bc5SJunchao Zhang 270d5b485abSSatish Balay /* receive work done on other processors*/ 271d5b485abSSatish Balay { 272b24ad042SBarry Smith PetscMPIInt idex; 273b24ad042SBarry Smith PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 274f1af5d2fSBarry Smith PetscBT table_i; 275d5b485abSSatish Balay 276d5b485abSSatish Balay for (i=0; i<nrqs; ++i) { 277175c2bc5SJunchao Zhang ierr = MPI_Waitany(nrqs,r_waits2,&idex,MPI_STATUS_IGNORE);CHKERRMPI(ierr); 278d5b485abSSatish Balay /* Process the message*/ 279999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 280d5b485abSSatish Balay ct1 = 2*rbuf2_i[0]+1; 281999d9058SBarry Smith jmax = rbuf2[idex][0]; 282d5b485abSSatish Balay for (j=1; j<=jmax; j++) { 283d5b485abSSatish Balay max = rbuf2_i[2*j]; 284d5b485abSSatish Balay is_no = rbuf2_i[2*j-1]; 285d5b485abSSatish Balay isz_i = isz[is_no]; 286d5b485abSSatish Balay data_i = data[is_no]; 287d5b485abSSatish Balay table_i = table[is_no]; 288d5b485abSSatish Balay for (k=0; k<max; k++,ct1++) { 289d5b485abSSatish Balay row = rbuf2_i[ct1]; 29026fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 291d5b485abSSatish Balay } 292d5b485abSSatish Balay isz[is_no] = isz_i; 293d5b485abSSatish Balay } 294d5b485abSSatish Balay } 295175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqr,s_waits2,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 296d5b485abSSatish Balay } 297d5b485abSSatish Balay 298d5b485abSSatish Balay for (i=0; i<imax; ++i) { 299749130a8SStefano Zampini ierr = ISCreateGeneral(iscomms[i],isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 300749130a8SStefano Zampini ierr = PetscCommDestroy(&iscomms[i]);CHKERRQ(ierr); 301d5b485abSSatish Balay } 302d5b485abSSatish Balay 303749130a8SStefano Zampini ierr = PetscFree(iscomms);CHKERRQ(ierr); 304c7dd2924SSatish Balay ierr = PetscFree(onodes2);CHKERRQ(ierr); 305c7dd2924SSatish Balay ierr = PetscFree(olengths2);CHKERRQ(ierr); 306c7dd2924SSatish Balay 307606d414cSSatish Balay ierr = PetscFree(pa);CHKERRQ(ierr); 308175c2bc5SJunchao Zhang if (rbuf2) { 309c05d87d6SBarry Smith ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 310606d414cSSatish Balay ierr = PetscFree(rbuf2);CHKERRQ(ierr); 311175c2bc5SJunchao Zhang } 312606d414cSSatish Balay ierr = PetscFree(s_waits1);CHKERRQ(ierr); 313606d414cSSatish Balay ierr = PetscFree(r_waits1);CHKERRQ(ierr); 314606d414cSSatish Balay ierr = PetscFree(s_waits2);CHKERRQ(ierr); 315606d414cSSatish Balay ierr = PetscFree(r_waits2);CHKERRQ(ierr); 3168f9f447aSHong Zhang ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 317175c2bc5SJunchao Zhang if (xdata) { 318606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 319606d414cSSatish Balay ierr = PetscFree(xdata);CHKERRQ(ierr); 320175c2bc5SJunchao Zhang } 321606d414cSSatish Balay ierr = PetscFree(isz1);CHKERRQ(ierr); 3223a40ed3dSBarry Smith PetscFunctionReturn(0); 323d5b485abSSatish Balay } 324d5b485abSSatish Balay 325d5b485abSSatish Balay /* 326df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 327d5b485abSSatish Balay the work on the local processor. 328d5b485abSSatish Balay 329d5b485abSSatish Balay Inputs: 330df36731bSBarry Smith C - MAT_MPIBAIJ; 331d5b485abSSatish Balay imax - total no of index sets processed at a time; 332df36731bSBarry Smith table - an array of char - size = Mbs bits. 333d5b485abSSatish Balay 334d5b485abSSatish Balay Output: 3350e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 336d5b485abSSatish Balay to each index set; 337d5b485abSSatish Balay data - pointer to the solutions 338d5b485abSSatish Balay */ 339b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 340d5b485abSSatish Balay { 341df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 342d5b485abSSatish Balay Mat A = c->A,B = c->B; 343df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 344b24ad042SBarry Smith PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 345b24ad042SBarry Smith PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 346f1af5d2fSBarry Smith PetscBT table_i; 347d5b485abSSatish Balay 3483a40ed3dSBarry Smith PetscFunctionBegin; 349899cda47SBarry Smith rstart = c->rstartbs; 350899cda47SBarry Smith cstart = c->cstartbs; 351d5b485abSSatish Balay ai = a->i; 352df36731bSBarry Smith aj = a->j; 353d5b485abSSatish Balay bi = b->i; 354df36731bSBarry Smith bj = b->j; 355d5b485abSSatish Balay garray = c->garray; 356d5b485abSSatish Balay 357d5b485abSSatish Balay for (i=0; i<imax; i++) { 358d5b485abSSatish Balay data_i = data[i]; 359d5b485abSSatish Balay table_i = table[i]; 360d5b485abSSatish Balay isz_i = isz[i]; 361d5b485abSSatish Balay for (j=0,max=isz[i]; j<max; j++) { 362d5b485abSSatish Balay row = data_i[j] - rstart; 363d5b485abSSatish Balay start = ai[row]; 364d5b485abSSatish Balay end = ai[row+1]; 365d5b485abSSatish Balay for (k=start; k<end; k++) { /* Amat */ 366df36731bSBarry Smith val = aj[k] + cstart; 36726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 368d5b485abSSatish Balay } 369d5b485abSSatish Balay start = bi[row]; 370d5b485abSSatish Balay end = bi[row+1]; 371d5b485abSSatish Balay for (k=start; k<end; k++) { /* Bmat */ 372df36731bSBarry Smith val = garray[bj[k]]; 37326fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 374d5b485abSSatish Balay } 375d5b485abSSatish Balay } 376d5b485abSSatish Balay isz[i] = isz_i; 377d5b485abSSatish Balay } 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 379d5b485abSSatish Balay } 380d5b485abSSatish Balay /* 3812d4ee042Sprj- MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages, 382d5b485abSSatish Balay and return the output 383d5b485abSSatish Balay 384d5b485abSSatish Balay Input: 385d5b485abSSatish Balay C - the matrix 386d5b485abSSatish Balay nrqr - no of messages being processed. 3872d4ee042Sprj- rbuf - an array of pointers to the received requests 388d5b485abSSatish Balay 389d5b485abSSatish Balay Output: 390d5b485abSSatish Balay xdata - array of messages to be sent back 391d5b485abSSatish Balay isz1 - size of each message 392d5b485abSSatish Balay 393a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 394d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 395a5b23f4aSJose E. Roman rather than all previous rows as it is now where a single large chunk of 396d5b485abSSatish Balay memory is used. 397d5b485abSSatish Balay 398d5b485abSSatish Balay */ 399b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 400d5b485abSSatish Balay { 401df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 402d5b485abSSatish Balay Mat A = c->A,B = c->B; 403df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 4046849ba73SBarry Smith PetscErrorCode ierr; 405b24ad042SBarry Smith PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 406b24ad042SBarry Smith PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 407d2a212eaSBarry Smith PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 408b24ad042SBarry Smith PetscInt *rbuf_i,kmax,rbuf_0; 409f1af5d2fSBarry Smith PetscBT xtable; 410d5b485abSSatish Balay 4113a40ed3dSBarry Smith PetscFunctionBegin; 412df36731bSBarry Smith Mbs = c->Mbs; 413899cda47SBarry Smith rstart = c->rstartbs; 414899cda47SBarry Smith cstart = c->cstartbs; 415d5b485abSSatish Balay ai = a->i; 416df36731bSBarry Smith aj = a->j; 417d5b485abSSatish Balay bi = b->i; 418df36731bSBarry Smith bj = b->j; 419d5b485abSSatish Balay garray = c->garray; 420d5b485abSSatish Balay 421d5b485abSSatish Balay for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 422d5b485abSSatish Balay rbuf_i = rbuf[i]; 423d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 424d5b485abSSatish Balay ct += rbuf_0; 42526fbe8dcSKarl Rupp for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 426d5b485abSSatish Balay } 427d5b485abSSatish Balay 428701b8814SBarry Smith if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 429701b8814SBarry Smith else max1 = 1; 430d5b485abSSatish Balay mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 431175c2bc5SJunchao Zhang if (nrqr) { 432785e854fSJed Brown ierr = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr); 433d5b485abSSatish Balay ++no_malloc; 434175c2bc5SJunchao Zhang } 43553b8de81SBarry Smith ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 436580bdb30SBarry Smith ierr = PetscArrayzero(isz1,nrqr);CHKERRQ(ierr); 437d5b485abSSatish Balay 438d5b485abSSatish Balay ct3 = 0; 439d5b485abSSatish Balay for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 440d5b485abSSatish Balay rbuf_i = rbuf[i]; 441d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 442d5b485abSSatish Balay ct1 = 2*rbuf_0+1; 443d5b485abSSatish Balay ct2 = ct1; 444d5b485abSSatish Balay ct3 += ct1; 445d5b485abSSatish Balay for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 4466831982aSBarry Smith ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 447d5b485abSSatish Balay oct2 = ct2; 448d5b485abSSatish Balay kmax = rbuf_i[2*j]; 449d5b485abSSatish Balay for (k=0; k<kmax; k++,ct1++) { 450d5b485abSSatish Balay row = rbuf_i[ct1]; 451f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,row)) { 452d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 453b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 454785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 455580bdb30SBarry Smith ierr = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr); 456606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 457d5b485abSSatish Balay xdata[0] = tmp; 458d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 45926fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 460d5b485abSSatish Balay } 461d5b485abSSatish Balay xdata[i][ct2++] = row; 462d5b485abSSatish Balay ct3++; 463d5b485abSSatish Balay } 464d5b485abSSatish Balay } 465d5b485abSSatish Balay for (k=oct2,max2=ct2; k<max2; k++) { 466d5b485abSSatish Balay row = xdata[i][k] - rstart; 467d5b485abSSatish Balay start = ai[row]; 468d5b485abSSatish Balay end = ai[row+1]; 469d5b485abSSatish Balay for (l=start; l<end; l++) { 470df36731bSBarry Smith val = aj[l] + cstart; 471f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 472d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 473b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 474785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 475580bdb30SBarry Smith ierr = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr); 476606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 477d5b485abSSatish Balay xdata[0] = tmp; 478d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 47926fbe8dcSKarl Rupp for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 480d5b485abSSatish Balay } 481d5b485abSSatish Balay xdata[i][ct2++] = val; 482d5b485abSSatish Balay ct3++; 483d5b485abSSatish Balay } 484d5b485abSSatish Balay } 485d5b485abSSatish Balay start = bi[row]; 486d5b485abSSatish Balay end = bi[row+1]; 487d5b485abSSatish Balay for (l=start; l<end; l++) { 488df36731bSBarry Smith val = garray[bj[l]]; 489f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable,val)) { 490d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 491b24ad042SBarry Smith new_estimate = (PetscInt)(1.5*mem_estimate)+1; 492785e854fSJed Brown ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 493580bdb30SBarry Smith ierr = PetscArraycpy(tmp,xdata[0],mem_estimate);CHKERRQ(ierr); 494606d414cSSatish Balay ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 495d5b485abSSatish Balay xdata[0] = tmp; 496d5b485abSSatish Balay mem_estimate = new_estimate; ++no_malloc; 49726fbe8dcSKarl Rupp for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 498d5b485abSSatish Balay } 499d5b485abSSatish Balay xdata[i][ct2++] = val; 500d5b485abSSatish Balay ct3++; 501d5b485abSSatish Balay } 502d5b485abSSatish Balay } 503d5b485abSSatish Balay } 504d5b485abSSatish Balay /* Update the header*/ 505d5b485abSSatish Balay xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 506d5b485abSSatish Balay xdata[i][2*j-1] = rbuf_i[2*j-1]; 507d5b485abSSatish Balay } 508d5b485abSSatish Balay xdata[i][0] = rbuf_0; 509175c2bc5SJunchao Zhang if (i+1<nrqr) xdata[i+1] = xdata[i] + ct2; 510d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 511d5b485abSSatish Balay } 51294bacf5dSBarry Smith ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 513*7d3de750SJacob Faibussowitsch ierr = PetscInfo(C,"Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 5143a40ed3dSBarry Smith PetscFunctionReturn(0); 515d5b485abSSatish Balay } 516d5b485abSSatish Balay 5177dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 518d5b485abSSatish Balay { 51917df9f7cSHong Zhang IS *isrow_block,*iscol_block; 520cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 5216849ba73SBarry Smith PetscErrorCode ierr; 522121971b2SHong Zhang PetscInt nmax,nstages_local,nstages,i,pos,max_no,N=C->cmap->N,bs=C->rmap->bs; 5235dc98d6aSHong Zhang Mat_SeqBAIJ *subc; 5245c39f6d9SHong Zhang Mat_SubSppt *smat; 525a2feddc5SSatish Balay 5263a40ed3dSBarry Smith PetscFunctionBegin; 52729dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 52829dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 52917df9f7cSHong Zhang ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr); 53017df9f7cSHong Zhang ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr); 53117df9f7cSHong Zhang ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr); 532cf4f527aSSatish Balay 533cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 5344698041cSHong Zhang if (!C->cmap->N) nmax=20*1000000/sizeof(PetscInt); 5354698041cSHong Zhang else nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 536329f5518SBarry Smith if (!nmax) nmax = 1; 5375dc98d6aSHong Zhang 5385dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 5394698041cSHong Zhang nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */ 540cf4f527aSSatish Balay 541653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 542820f2d46SBarry Smith ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRMPI(ierr); 5435dc98d6aSHong Zhang 5444698041cSHong Zhang /* Allocate memory to hold all the submatrices and dummy submatrices */ 5454698041cSHong Zhang ierr = PetscCalloc1(ismax+nstages,submat);CHKERRQ(ierr); 5464698041cSHong Zhang } else { /* MAT_REUSE_MATRIX */ 5474698041cSHong Zhang if (ismax) { 5485dc98d6aSHong Zhang subc = (Mat_SeqBAIJ*)((*submat)[0]->data); 5495dc98d6aSHong Zhang smat = subc->submatis1; 5504698041cSHong Zhang } else { /* (*submat)[0] is a dummy matrix */ 5515c39f6d9SHong Zhang smat = (Mat_SubSppt*)(*submat)[0]->data; 5524698041cSHong Zhang } 553071fcb05SBarry Smith if (!smat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat"); 5545dc98d6aSHong Zhang nstages = smat->nstages; 5555dc98d6aSHong Zhang } 5565dc98d6aSHong Zhang 557cf4f527aSSatish Balay for (i=0,pos=0; i<nstages; i++) { 558cf4f527aSSatish Balay if (pos+nmax <= ismax) max_no = nmax; 5594e195584SJunchao Zhang else if (pos >= ismax) max_no = 0; 560cf4f527aSSatish Balay else max_no = ismax-pos; 561a42ab0d6SHong Zhang 5627dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos);CHKERRQ(ierr); 5634e195584SJunchao Zhang if (!max_no) { 5644e195584SJunchao Zhang if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */ 5655c39f6d9SHong Zhang smat = (Mat_SubSppt*)(*submat)[pos]->data; 5664698041cSHong Zhang smat->nstages = nstages; 5674698041cSHong Zhang } 5684e195584SJunchao Zhang pos++; /* advance to next dummy matrix if any */ 5694e195584SJunchao Zhang } else pos += max_no; 570cf4f527aSSatish Balay } 57136f4e84dSSatish Balay 5725dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX && ismax) { 5735dc98d6aSHong Zhang /* save nstages for reuse */ 5745dc98d6aSHong Zhang subc = (Mat_SeqBAIJ*)((*submat)[0]->data); 5755dc98d6aSHong Zhang smat = subc->submatis1; 5765dc98d6aSHong Zhang smat->nstages = nstages; 5775dc98d6aSHong Zhang } 5785dc98d6aSHong Zhang 57936f4e84dSSatish Balay for (i=0; i<ismax; i++) { 58017df9f7cSHong Zhang ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr); 58117df9f7cSHong Zhang ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr); 58236f4e84dSSatish Balay } 58317df9f7cSHong Zhang ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr); 5843a40ed3dSBarry Smith PetscFunctionReturn(0); 585a2feddc5SSatish Balay } 586a2feddc5SSatish Balay 587233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 588e005ede5SBarry Smith PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 589233c2ff6SSatish Balay { 590e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 5914dc2109aSBarry Smith PetscMPIInt fproc; 5924dc2109aSBarry Smith PetscErrorCode ierr; 593233c2ff6SSatish Balay 594233c2ff6SSatish Balay PetscFunctionBegin; 5954dc2109aSBarry Smith ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 59623ce1328SBarry Smith if (fproc > size) fproc = size; 597e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 598e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 599233c2ff6SSatish Balay else fproc++; 600233c2ff6SSatish Balay } 601e005ede5SBarry Smith *rank = fproc; 602233c2ff6SSatish Balay PetscFunctionReturn(0); 603233c2ff6SSatish Balay } 604233c2ff6SSatish Balay #endif 605233c2ff6SSatish Balay 606a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 607b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6087dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 60917df9f7cSHong Zhang { 61017df9f7cSHong Zhang Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 61117df9f7cSHong Zhang Mat A = c->A; 61217df9f7cSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc; 61317df9f7cSHong Zhang const PetscInt **icol,**irow; 61417df9f7cSHong Zhang PetscInt *nrow,*ncol,start; 61517df9f7cSHong Zhang PetscErrorCode ierr; 61617df9f7cSHong Zhang PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr; 617ddb3d6bcSHong Zhang PetscInt **sbuf1,**sbuf2,*sbuf2_i,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1; 61817df9f7cSHong Zhang PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol; 61917df9f7cSHong Zhang PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2; 62017df9f7cSHong Zhang PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 62117df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 62217df9f7cSHong Zhang PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i; 62317df9f7cSHong Zhang #else 62417df9f7cSHong Zhang PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i; 62517df9f7cSHong Zhang #endif 62696cff833SHong Zhang const PetscInt *irow_i,*icol_i; 62717df9f7cSHong Zhang PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 62817df9f7cSHong Zhang MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 62917df9f7cSHong Zhang MPI_Request *r_waits4,*s_waits3,*s_waits4; 63017df9f7cSHong Zhang MPI_Comm comm; 63164f5bb2dSSatish Balay PetscScalar **rbuf4,*rbuf4_i=NULL,**sbuf_aa,*vals,*mat_a=NULL,*imat_a=NULL,*sbuf_aa_i; 63217df9f7cSHong Zhang PetscMPIInt *onodes1,*olengths1,end; 63380d31651SHong Zhang PetscInt **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i; 6345c39f6d9SHong Zhang Mat_SubSppt *smat_i; 63517df9f7cSHong Zhang PetscBool *issorted,colflag,iscsorted=PETSC_TRUE; 63617df9f7cSHong Zhang PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen; 637a42ab0d6SHong Zhang PetscInt bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs; 638a42ab0d6SHong Zhang PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 639a42ab0d6SHong Zhang PetscInt nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB; 640afb3cbd8SHong Zhang PetscScalar *vworkA=NULL,*vworkB=NULL,*a_a = a->a,*b_a = b->a; 641a42ab0d6SHong Zhang PetscInt cstart = c->cstartbs,*bmap = c->garray; 64280d31651SHong Zhang PetscBool *allrows,*allcolumns; 643a42ab0d6SHong Zhang 64417df9f7cSHong Zhang PetscFunctionBegin; 64517df9f7cSHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 64617df9f7cSHong Zhang size = c->size; 64717df9f7cSHong Zhang rank = c->rank; 64817df9f7cSHong Zhang 6494698041cSHong Zhang ierr = PetscMalloc5(ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax+1,&allcolumns,ismax,&allrows);CHKERRQ(ierr); 650f489ac74SBarry Smith ierr = PetscMalloc5(ismax,(PetscInt***)&irow,ismax,(PetscInt***)&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr); 65117df9f7cSHong Zhang 65217df9f7cSHong Zhang for (i=0; i<ismax; i++) { 65317df9f7cSHong Zhang ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr); 65480d31651SHong Zhang if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */ 65517df9f7cSHong Zhang ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr); 65617df9f7cSHong Zhang 657ec3c739cSHong Zhang /* Check for special case: allcolumns */ 65817df9f7cSHong Zhang ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 65917df9f7cSHong Zhang ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 6605dd0c08cSHong Zhang 6611abc2f7bSHong Zhang if (colflag && ncol[i] == c->Nbs) { 6621abc2f7bSHong Zhang allcolumns[i] = PETSC_TRUE; 6631abc2f7bSHong Zhang icol[i] = NULL; 6641abc2f7bSHong Zhang } else { 66517df9f7cSHong Zhang allcolumns[i] = PETSC_FALSE; 66617df9f7cSHong Zhang ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 6671abc2f7bSHong Zhang } 668ec3c739cSHong Zhang 669ec3c739cSHong Zhang /* Check for special case: allrows */ 670ec3c739cSHong Zhang ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 671ec3c739cSHong Zhang ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 672ec3c739cSHong Zhang if (colflag && nrow[i] == c->Mbs) { 673ec3c739cSHong Zhang allrows[i] = PETSC_TRUE; 674ec3c739cSHong Zhang irow[i] = NULL; 675ec3c739cSHong Zhang } else { 676ec3c739cSHong Zhang allrows[i] = PETSC_FALSE; 677ec3c739cSHong Zhang ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 678ec3c739cSHong Zhang } 67917df9f7cSHong Zhang } 68017df9f7cSHong Zhang 68117df9f7cSHong Zhang if (scall == MAT_REUSE_MATRIX) { 68217df9f7cSHong Zhang /* Assumes new rows are same length as the old rows */ 68317df9f7cSHong Zhang for (i=0; i<ismax; i++) { 68417df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)(submats[i]->data); 6859e686edcSHong Zhang if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 68617df9f7cSHong Zhang 68717df9f7cSHong Zhang /* Initial matrix as if empty */ 688580bdb30SBarry Smith ierr = PetscArrayzero(subc->ilen,subc->mbs);CHKERRQ(ierr); 68917df9f7cSHong Zhang 69017df9f7cSHong Zhang /* Initial matrix as if empty */ 69117df9f7cSHong Zhang submats[i]->factortype = C->factortype; 69217df9f7cSHong Zhang 69317df9f7cSHong Zhang smat_i = subc->submatis1; 69417df9f7cSHong Zhang 69517df9f7cSHong Zhang nrqs = smat_i->nrqs; 69617df9f7cSHong Zhang nrqr = smat_i->nrqr; 69717df9f7cSHong Zhang rbuf1 = smat_i->rbuf1; 69817df9f7cSHong Zhang rbuf2 = smat_i->rbuf2; 69917df9f7cSHong Zhang rbuf3 = smat_i->rbuf3; 70017df9f7cSHong Zhang req_source2 = smat_i->req_source2; 70117df9f7cSHong Zhang 70217df9f7cSHong Zhang sbuf1 = smat_i->sbuf1; 70317df9f7cSHong Zhang sbuf2 = smat_i->sbuf2; 70417df9f7cSHong Zhang ptr = smat_i->ptr; 70517df9f7cSHong Zhang tmp = smat_i->tmp; 70617df9f7cSHong Zhang ctr = smat_i->ctr; 70717df9f7cSHong Zhang 70817df9f7cSHong Zhang pa = smat_i->pa; 70917df9f7cSHong Zhang req_size = smat_i->req_size; 71017df9f7cSHong Zhang req_source1 = smat_i->req_source1; 71117df9f7cSHong Zhang 71217df9f7cSHong Zhang allcolumns[i] = smat_i->allcolumns; 713ec3c739cSHong Zhang allrows[i] = smat_i->allrows; 71417df9f7cSHong Zhang row2proc[i] = smat_i->row2proc; 71517df9f7cSHong Zhang rmap[i] = smat_i->rmap; 71617df9f7cSHong Zhang cmap[i] = smat_i->cmap; 71717df9f7cSHong Zhang } 7184698041cSHong Zhang 7194698041cSHong Zhang if (!ismax){ /* Get dummy submatrices and retrieve struct submatis1 */ 7204698041cSHong Zhang if (!submats[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats are null, cannot reuse"); 7215c39f6d9SHong Zhang smat_i = (Mat_SubSppt*)submats[0]->data; 7224698041cSHong Zhang 7234698041cSHong Zhang nrqs = smat_i->nrqs; 7244698041cSHong Zhang nrqr = smat_i->nrqr; 7254698041cSHong Zhang rbuf1 = smat_i->rbuf1; 7264698041cSHong Zhang rbuf2 = smat_i->rbuf2; 7274698041cSHong Zhang rbuf3 = smat_i->rbuf3; 7284698041cSHong Zhang req_source2 = smat_i->req_source2; 7294698041cSHong Zhang 7304698041cSHong Zhang sbuf1 = smat_i->sbuf1; 7314698041cSHong Zhang sbuf2 = smat_i->sbuf2; 7324698041cSHong Zhang ptr = smat_i->ptr; 7334698041cSHong Zhang tmp = smat_i->tmp; 7344698041cSHong Zhang ctr = smat_i->ctr; 7354698041cSHong Zhang 7364698041cSHong Zhang pa = smat_i->pa; 7374698041cSHong Zhang req_size = smat_i->req_size; 7384698041cSHong Zhang req_source1 = smat_i->req_source1; 7394698041cSHong Zhang 740c67fe082SHong Zhang allcolumns[0] = PETSC_FALSE; 7414698041cSHong Zhang } 74217df9f7cSHong Zhang } else { /* scall == MAT_INITIAL_MATRIX */ 74317df9f7cSHong Zhang /* Get some new tags to keep the communication clean */ 74417df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 74517df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 74617df9f7cSHong Zhang 74717df9f7cSHong Zhang /* evaluate communication - mesg to who, length of mesg, and buffer space 74817df9f7cSHong Zhang required. Based on this, buffers are allocated, and data copied into them*/ 74917df9f7cSHong Zhang ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); /* mesg size, initialize work vectors */ 75017df9f7cSHong Zhang 75117df9f7cSHong Zhang for (i=0; i<ismax; i++) { 75217df9f7cSHong Zhang jmax = nrow[i]; 75317df9f7cSHong Zhang irow_i = irow[i]; 75417df9f7cSHong Zhang 75517df9f7cSHong Zhang ierr = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr); 75617df9f7cSHong Zhang row2proc[i] = row2proc_i; 75717df9f7cSHong Zhang 75817df9f7cSHong Zhang if (issorted[i]) proc = 0; 75917df9f7cSHong Zhang for (j=0; j<jmax; j++) { 76017df9f7cSHong Zhang if (!issorted[i]) proc = 0; 761ec3c739cSHong Zhang if (allrows[i]) row = j; 762ec3c739cSHong Zhang else row = irow_i[j]; 763ec3c739cSHong Zhang 764a42ab0d6SHong Zhang while (row >= c->rangebs[proc+1]) proc++; 76517df9f7cSHong Zhang w4[proc]++; 76617df9f7cSHong Zhang row2proc_i[j] = proc; /* map row index to proc */ 76717df9f7cSHong Zhang } 76817df9f7cSHong Zhang for (j=0; j<size; j++) { 76917df9f7cSHong Zhang if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;} 77017df9f7cSHong Zhang } 77117df9f7cSHong Zhang } 77217df9f7cSHong Zhang 77317df9f7cSHong Zhang nrqs = 0; /* no of outgoing messages */ 77417df9f7cSHong Zhang msz = 0; /* total mesg length (for all procs) */ 77517df9f7cSHong Zhang w1[rank] = 0; /* no mesg sent to self */ 77617df9f7cSHong Zhang w3[rank] = 0; 77717df9f7cSHong Zhang for (i=0; i<size; i++) { 77817df9f7cSHong Zhang if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 77917df9f7cSHong Zhang } 780175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&pa);CHKERRQ(ierr); /*(proc -array)*/ 78117df9f7cSHong Zhang for (i=0,j=0; i<size; i++) { 78217df9f7cSHong Zhang if (w1[i]) { pa[j] = i; j++; } 78317df9f7cSHong Zhang } 78417df9f7cSHong Zhang 78517df9f7cSHong Zhang /* Each message would have a header = 1 + 2*(no of IS) + data */ 78617df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 78717df9f7cSHong Zhang j = pa[i]; 78817df9f7cSHong Zhang w1[j] += w2[j] + 2* w3[j]; 78917df9f7cSHong Zhang msz += w1[j]; 79017df9f7cSHong Zhang } 791*7d3de750SJacob Faibussowitsch ierr = PetscInfo(0,"Number of outgoing messages %" PetscInt_FMT " Total message length %" PetscInt_FMT "\n",nrqs,msz);CHKERRQ(ierr); 79217df9f7cSHong Zhang 79317df9f7cSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 79417df9f7cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 79517df9f7cSHong Zhang ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 79617df9f7cSHong Zhang 79717df9f7cSHong Zhang /* Now post the Irecvs corresponding to these messages */ 798175c2bc5SJunchao Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag0);CHKERRQ(ierr); 79917df9f7cSHong Zhang ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 80017df9f7cSHong Zhang 80117df9f7cSHong Zhang /* Allocate Memory for outgoing messages */ 80217df9f7cSHong Zhang ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 803580bdb30SBarry Smith ierr = PetscArrayzero(sbuf1,size);CHKERRQ(ierr); 804580bdb30SBarry Smith ierr = PetscArrayzero(ptr,size);CHKERRQ(ierr); 80517df9f7cSHong Zhang 80617df9f7cSHong Zhang { 80717df9f7cSHong Zhang PetscInt *iptr = tmp; 80817df9f7cSHong Zhang k = 0; 80917df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 81017df9f7cSHong Zhang j = pa[i]; 81117df9f7cSHong Zhang iptr += k; 81217df9f7cSHong Zhang sbuf1[j] = iptr; 81317df9f7cSHong Zhang k = w1[j]; 81417df9f7cSHong Zhang } 81517df9f7cSHong Zhang } 81617df9f7cSHong Zhang 81717df9f7cSHong Zhang /* Form the outgoing messages. Initialize the header space */ 81817df9f7cSHong Zhang for (i=0; i<nrqs; i++) { 81917df9f7cSHong Zhang j = pa[i]; 82017df9f7cSHong Zhang sbuf1[j][0] = 0; 821580bdb30SBarry Smith ierr = PetscArrayzero(sbuf1[j]+1,2*w3[j]);CHKERRQ(ierr); 82217df9f7cSHong Zhang ptr[j] = sbuf1[j] + 2*w3[j] + 1; 82317df9f7cSHong Zhang } 82417df9f7cSHong Zhang 82517df9f7cSHong Zhang /* Parse the isrow and copy data into outbuf */ 82617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 82717df9f7cSHong Zhang row2proc_i = row2proc[i]; 828580bdb30SBarry Smith ierr = PetscArrayzero(ctr,size);CHKERRQ(ierr); 82917df9f7cSHong Zhang irow_i = irow[i]; 83017df9f7cSHong Zhang jmax = nrow[i]; 83117df9f7cSHong Zhang for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 83217df9f7cSHong Zhang proc = row2proc_i[j]; 833ec3c739cSHong Zhang if (allrows[i]) row = j; 834ec3c739cSHong Zhang else row = irow_i[j]; 835ec3c739cSHong Zhang 83617df9f7cSHong Zhang if (proc != rank) { /* copy to the outgoing buf*/ 83717df9f7cSHong Zhang ctr[proc]++; 838ec3c739cSHong Zhang *ptr[proc] = row; 83917df9f7cSHong Zhang ptr[proc]++; 84017df9f7cSHong Zhang } 84117df9f7cSHong Zhang } 84217df9f7cSHong Zhang /* Update the headers for the current IS */ 84317df9f7cSHong Zhang for (j=0; j<size; j++) { /* Can Optimise this loop too */ 84417df9f7cSHong Zhang if ((ctr_j = ctr[j])) { 84517df9f7cSHong Zhang sbuf1_j = sbuf1[j]; 84617df9f7cSHong Zhang k = ++sbuf1_j[0]; 84717df9f7cSHong Zhang sbuf1_j[2*k] = ctr_j; 84817df9f7cSHong Zhang sbuf1_j[2*k-1] = i; 84917df9f7cSHong Zhang } 85017df9f7cSHong Zhang } 85117df9f7cSHong Zhang } 85217df9f7cSHong Zhang 85317df9f7cSHong Zhang /* Now post the sends */ 854175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&s_waits1);CHKERRQ(ierr); 85517df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 85617df9f7cSHong Zhang j = pa[i]; 857ffc4695bSBarry Smith ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRMPI(ierr); 85817df9f7cSHong Zhang } 85917df9f7cSHong Zhang 86017df9f7cSHong Zhang /* Post Receives to capture the buffer size */ 861175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&r_waits2);CHKERRQ(ierr); 862175c2bc5SJunchao Zhang ierr = PetscMalloc3(nrqs,&req_source2,nrqs,&rbuf2,nrqs,&rbuf3);CHKERRQ(ierr); 863175c2bc5SJunchao Zhang if (nrqs) rbuf2[0] = tmp + msz; 86417df9f7cSHong Zhang for (i=1; i<nrqs; ++i) { 86517df9f7cSHong Zhang rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 86617df9f7cSHong Zhang } 86717df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 86817df9f7cSHong Zhang j = pa[i]; 869ffc4695bSBarry Smith ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRMPI(ierr); 87017df9f7cSHong Zhang } 87117df9f7cSHong Zhang 87217df9f7cSHong Zhang /* Send to other procs the buf size they should allocate */ 87317df9f7cSHong Zhang /* Receive messages*/ 874175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&s_waits2);CHKERRQ(ierr); 87517df9f7cSHong Zhang ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr); 87617df9f7cSHong Zhang 877175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqr,r_waits1,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 87817df9f7cSHong Zhang for (i=0; i<nrqr; ++i) { 87917df9f7cSHong Zhang req_size[i] = 0; 88017df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 88117df9f7cSHong Zhang start = 2*rbuf1_i[0] + 1; 882175c2bc5SJunchao Zhang end = olengths1[i]; 883175c2bc5SJunchao Zhang ierr = PetscMalloc1(end,&sbuf2[i]);CHKERRQ(ierr); 88417df9f7cSHong Zhang sbuf2_i = sbuf2[i]; 88517df9f7cSHong Zhang for (j=start; j<end; j++) { 886ddb3d6bcSHong Zhang row = rbuf1_i[j] - rstart; 887ddb3d6bcSHong Zhang ncols = a_i[row+1] - a_i[row] + b_i[row+1] - b_i[row]; 88817df9f7cSHong Zhang sbuf2_i[j] = ncols; 88917df9f7cSHong Zhang req_size[i] += ncols; 89017df9f7cSHong Zhang } 891175c2bc5SJunchao Zhang req_source1[i] = onodes1[i]; 89217df9f7cSHong Zhang /* form the header */ 89317df9f7cSHong Zhang sbuf2_i[0] = req_size[i]; 89417df9f7cSHong Zhang for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 89517df9f7cSHong Zhang 896ffc4695bSBarry Smith ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRMPI(ierr); 89717df9f7cSHong Zhang } 898ddb3d6bcSHong Zhang 899175c2bc5SJunchao Zhang ierr = PetscFree(onodes1);CHKERRQ(ierr); 900175c2bc5SJunchao Zhang ierr = PetscFree(olengths1);CHKERRQ(ierr); 901175c2bc5SJunchao Zhang 90217df9f7cSHong Zhang ierr = PetscFree(r_waits1);CHKERRQ(ierr); 90317df9f7cSHong Zhang ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 90417df9f7cSHong Zhang 90517df9f7cSHong Zhang /* Receive messages*/ 906175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&r_waits3);CHKERRQ(ierr); 90717df9f7cSHong Zhang 908175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqs,r_waits2,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 90917df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 910175c2bc5SJunchao Zhang ierr = PetscMalloc1(rbuf2[i][0],&rbuf3[i]);CHKERRQ(ierr); 911175c2bc5SJunchao Zhang req_source2[i] = pa[i]; 912ffc4695bSBarry Smith ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRMPI(ierr); 91317df9f7cSHong Zhang } 91417df9f7cSHong Zhang ierr = PetscFree(r_waits2);CHKERRQ(ierr); 91517df9f7cSHong Zhang 91617df9f7cSHong Zhang /* Wait on sends1 and sends2 */ 917175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqs,s_waits1,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 918175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqr,s_waits2,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 91917df9f7cSHong Zhang ierr = PetscFree(s_waits1);CHKERRQ(ierr); 92017df9f7cSHong Zhang ierr = PetscFree(s_waits2);CHKERRQ(ierr); 92117df9f7cSHong Zhang 92217df9f7cSHong Zhang /* Now allocate sending buffers for a->j, and send them off */ 923175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&sbuf_aj);CHKERRQ(ierr); 92417df9f7cSHong Zhang for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 925175c2bc5SJunchao Zhang if (nrqr) {ierr = PetscMalloc1(j,&sbuf_aj[0]);CHKERRQ(ierr);} 92617df9f7cSHong Zhang for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 92717df9f7cSHong Zhang 928175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&s_waits3);CHKERRQ(ierr); 92917df9f7cSHong Zhang { 93017df9f7cSHong Zhang 93117df9f7cSHong Zhang for (i=0; i<nrqr; i++) { 93217df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 93317df9f7cSHong Zhang sbuf_aj_i = sbuf_aj[i]; 93417df9f7cSHong Zhang ct1 = 2*rbuf1_i[0] + 1; 93517df9f7cSHong Zhang ct2 = 0; 93617df9f7cSHong Zhang for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 93717df9f7cSHong Zhang kmax = rbuf1[i][2*j]; 93817df9f7cSHong Zhang for (k=0; k<kmax; k++,ct1++) { 93917df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 94017df9f7cSHong Zhang nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 94117df9f7cSHong Zhang ncols = nzA + nzB; 94217df9f7cSHong Zhang cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 94317df9f7cSHong Zhang 94417df9f7cSHong Zhang /* load the column indices for this row into cols */ 94517df9f7cSHong Zhang cols = sbuf_aj_i + ct2; 94617df9f7cSHong Zhang for (l=0; l<nzB; l++) { 947a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 948a42ab0d6SHong Zhang else break; 94917df9f7cSHong Zhang } 950a42ab0d6SHong Zhang imark = l; 951a42ab0d6SHong Zhang for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];} 952a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 95317df9f7cSHong Zhang ct2 += ncols; 95417df9f7cSHong Zhang } 95517df9f7cSHong Zhang } 956ffc4695bSBarry Smith ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRMPI(ierr); 95717df9f7cSHong Zhang } 95817df9f7cSHong Zhang } 95917df9f7cSHong Zhang 96017df9f7cSHong Zhang /* create col map: global col of C -> local col of submatrices */ 96117df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 96217df9f7cSHong Zhang for (i=0; i<ismax; i++) { 96317df9f7cSHong Zhang if (!allcolumns[i]) { 964175c2bc5SJunchao Zhang ierr = PetscTableCreate(ncol[i],c->Nbs,&cmap[i]);CHKERRQ(ierr); 96517df9f7cSHong Zhang 96617df9f7cSHong Zhang jmax = ncol[i]; 96717df9f7cSHong Zhang icol_i = icol[i]; 96817df9f7cSHong Zhang cmap_i = cmap[i]; 96917df9f7cSHong Zhang for (j=0; j<jmax; j++) { 97017df9f7cSHong Zhang ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 97117df9f7cSHong Zhang } 97217df9f7cSHong Zhang } else cmap[i] = NULL; 97317df9f7cSHong Zhang } 97417df9f7cSHong Zhang #else 97517df9f7cSHong Zhang for (i=0; i<ismax; i++) { 97617df9f7cSHong Zhang if (!allcolumns[i]) { 977a42ab0d6SHong Zhang ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 97817df9f7cSHong Zhang jmax = ncol[i]; 97917df9f7cSHong Zhang icol_i = icol[i]; 98017df9f7cSHong Zhang cmap_i = cmap[i]; 981a42ab0d6SHong Zhang for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1; 98217df9f7cSHong Zhang } else cmap[i] = NULL; 98317df9f7cSHong Zhang } 98417df9f7cSHong Zhang #endif 98517df9f7cSHong Zhang 98617df9f7cSHong Zhang /* Create lens which is required for MatCreate... */ 98717df9f7cSHong Zhang for (i=0,j=0; i<ismax; i++) j += nrow[i]; 98817df9f7cSHong Zhang ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr); 98917df9f7cSHong Zhang 99017df9f7cSHong Zhang if (ismax) { 99117df9f7cSHong Zhang ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr); 99217df9f7cSHong Zhang } 99317df9f7cSHong Zhang for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 99417df9f7cSHong Zhang 99517df9f7cSHong Zhang /* Update lens from local data */ 99617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 99717df9f7cSHong Zhang row2proc_i = row2proc[i]; 99817df9f7cSHong Zhang jmax = nrow[i]; 99917df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 100017df9f7cSHong Zhang irow_i = irow[i]; 100117df9f7cSHong Zhang lens_i = lens[i]; 100217df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1003ec3c739cSHong Zhang if (allrows[i]) row = j; 1004ec3c739cSHong Zhang else row = irow_i[j]; /* global blocked row of C */ 1005ec3c739cSHong Zhang 100617df9f7cSHong Zhang proc = row2proc_i[j]; 100717df9f7cSHong Zhang if (proc == rank) { 1008a42ab0d6SHong Zhang /* Get indices from matA and then from matB */ 100917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1010a42ab0d6SHong Zhang PetscInt tt; 1011a42ab0d6SHong Zhang #endif 1012a42ab0d6SHong Zhang row = row - rstart; 1013a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1014a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 1015a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1016a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1017a42ab0d6SHong Zhang 1018a42ab0d6SHong Zhang if (!allcolumns[i]) { 1019a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1020a42ab0d6SHong Zhang for (k=0; k<nzA; k++) { 1021a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1022a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1023a42ab0d6SHong Zhang } 1024a42ab0d6SHong Zhang for (k=0; k<nzB; k++) { 1025a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1026a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1027a42ab0d6SHong Zhang } 1028a42ab0d6SHong Zhang 1029a42ab0d6SHong Zhang #else 1030a42ab0d6SHong Zhang for (k=0; k<nzA; k++) { 1031a42ab0d6SHong Zhang if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1032a42ab0d6SHong Zhang } 1033a42ab0d6SHong Zhang for (k=0; k<nzB; k++) { 1034a42ab0d6SHong Zhang if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1035a42ab0d6SHong Zhang } 1036a42ab0d6SHong Zhang #endif 1037a42ab0d6SHong Zhang } else { /* allcolumns */ 1038a42ab0d6SHong Zhang lens_i[j] = nzA + nzB; 1039a42ab0d6SHong Zhang } 104017df9f7cSHong Zhang } 104117df9f7cSHong Zhang } 104217df9f7cSHong Zhang } 104317df9f7cSHong Zhang 104417df9f7cSHong Zhang /* Create row map: global row of C -> local row of submatrices */ 104517df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1046059f6b74SHong Zhang if (!allrows[i]) { 1047059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE) 1048175c2bc5SJunchao Zhang ierr = PetscTableCreate(nrow[i],c->Mbs,&rmap[i]);CHKERRQ(ierr); 104917df9f7cSHong Zhang irow_i = irow[i]; 105017df9f7cSHong Zhang jmax = nrow[i]; 105117df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1052ec3c739cSHong Zhang if (allrows[i]) { 1053ec3c739cSHong Zhang ierr = PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1054ec3c739cSHong Zhang } else { 105517df9f7cSHong Zhang ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 105617df9f7cSHong Zhang } 105717df9f7cSHong Zhang } 105817df9f7cSHong Zhang #else 1059a42ab0d6SHong Zhang ierr = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr); 106017df9f7cSHong Zhang rmap_i = rmap[i]; 106117df9f7cSHong Zhang irow_i = irow[i]; 106217df9f7cSHong Zhang jmax = nrow[i]; 106317df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1064ec3c739cSHong Zhang if (allrows[i]) rmap_i[j] = j; 1065ec3c739cSHong Zhang else rmap_i[irow_i[j]] = j; 106617df9f7cSHong Zhang } 106717df9f7cSHong Zhang #endif 1068059f6b74SHong Zhang } else rmap[i] = NULL; 1069059f6b74SHong Zhang } 107017df9f7cSHong Zhang 107117df9f7cSHong Zhang /* Update lens from offproc data */ 107217df9f7cSHong Zhang { 107317df9f7cSHong Zhang PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 107417df9f7cSHong Zhang 1075175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqs,r_waits3,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 107617df9f7cSHong Zhang for (tmp2=0; tmp2<nrqs; tmp2++) { 107717df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 107817df9f7cSHong Zhang jmax = sbuf1_i[0]; 107917df9f7cSHong Zhang ct1 = 2*jmax+1; 108017df9f7cSHong Zhang ct2 = 0; 108117df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 108217df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 108317df9f7cSHong Zhang for (j=1; j<=jmax; j++) { 108417df9f7cSHong Zhang is_no = sbuf1_i[2*j-1]; 108517df9f7cSHong Zhang max1 = sbuf1_i[2*j]; 108617df9f7cSHong Zhang lens_i = lens[is_no]; 108717df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 108817df9f7cSHong Zhang rmap_i = rmap[is_no]; 108917df9f7cSHong Zhang for (k=0; k<max1; k++,ct1++) { 1090059f6b74SHong Zhang if (allrows[is_no]) { 1091059f6b74SHong Zhang row = sbuf1_i[ct1]; 1092059f6b74SHong Zhang } else { 109317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 109417df9f7cSHong Zhang ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 109517df9f7cSHong Zhang row--; 109617df9f7cSHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 109717df9f7cSHong Zhang #else 109817df9f7cSHong Zhang row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 109917df9f7cSHong Zhang #endif 1100059f6b74SHong Zhang } 110117df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 110217df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 110317df9f7cSHong Zhang if (!allcolumns[is_no]) { 110417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 110517df9f7cSHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 110617df9f7cSHong Zhang #else 110717df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 110817df9f7cSHong Zhang #endif 110917df9f7cSHong Zhang if (tcol) lens_i[row]++; 111017df9f7cSHong Zhang } else { /* allcolumns */ 111117df9f7cSHong Zhang lens_i[row]++; /* lens_i[row] += max2 ? */ 111217df9f7cSHong Zhang } 111317df9f7cSHong Zhang } 111417df9f7cSHong Zhang } 111517df9f7cSHong Zhang } 111617df9f7cSHong Zhang } 111717df9f7cSHong Zhang } 111817df9f7cSHong Zhang ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1119175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqr,s_waits3,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 112017df9f7cSHong Zhang ierr = PetscFree(s_waits3);CHKERRQ(ierr); 112117df9f7cSHong Zhang 112217df9f7cSHong Zhang /* Create the submatrices */ 112317df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1124a42ab0d6SHong Zhang PetscInt bs_tmp; 1125a42ab0d6SHong Zhang if (ijonly) bs_tmp = 1; 1126a42ab0d6SHong Zhang else bs_tmp = bs; 112717df9f7cSHong Zhang 112817df9f7cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 1129a42ab0d6SHong Zhang ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 113017df9f7cSHong Zhang 113117df9f7cSHong Zhang ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 1132a42ab0d6SHong Zhang ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 1133a42ab0d6SHong Zhang ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 113417df9f7cSHong Zhang 11355c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 113617df9f7cSHong Zhang ierr = PetscNew(&smat_i);CHKERRQ(ierr); 113717df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 113817df9f7cSHong Zhang subc->submatis1 = smat_i; 113917df9f7cSHong Zhang 114017df9f7cSHong Zhang smat_i->destroy = submats[i]->ops->destroy; 1141f68bb481SHong Zhang submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ; 114217df9f7cSHong Zhang submats[i]->factortype = C->factortype; 114317df9f7cSHong Zhang 114417df9f7cSHong Zhang smat_i->id = i; 114517df9f7cSHong Zhang smat_i->nrqs = nrqs; 114617df9f7cSHong Zhang smat_i->nrqr = nrqr; 114717df9f7cSHong Zhang smat_i->rbuf1 = rbuf1; 114817df9f7cSHong Zhang smat_i->rbuf2 = rbuf2; 114917df9f7cSHong Zhang smat_i->rbuf3 = rbuf3; 115017df9f7cSHong Zhang smat_i->sbuf2 = sbuf2; 115117df9f7cSHong Zhang smat_i->req_source2 = req_source2; 115217df9f7cSHong Zhang 115317df9f7cSHong Zhang smat_i->sbuf1 = sbuf1; 115417df9f7cSHong Zhang smat_i->ptr = ptr; 115517df9f7cSHong Zhang smat_i->tmp = tmp; 115617df9f7cSHong Zhang smat_i->ctr = ctr; 115717df9f7cSHong Zhang 115817df9f7cSHong Zhang smat_i->pa = pa; 115917df9f7cSHong Zhang smat_i->req_size = req_size; 116017df9f7cSHong Zhang smat_i->req_source1 = req_source1; 116117df9f7cSHong Zhang 116217df9f7cSHong Zhang smat_i->allcolumns = allcolumns[i]; 1163ec3c739cSHong Zhang smat_i->allrows = allrows[i]; 116417df9f7cSHong Zhang smat_i->singleis = PETSC_FALSE; 116517df9f7cSHong Zhang smat_i->row2proc = row2proc[i]; 116617df9f7cSHong Zhang smat_i->rmap = rmap[i]; 116717df9f7cSHong Zhang smat_i->cmap = cmap[i]; 116817df9f7cSHong Zhang } 116917df9f7cSHong Zhang 11704698041cSHong Zhang if (!ismax) { /* Create dummy submats[0] for reuse struct subc */ 11714698041cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&submats[0]);CHKERRQ(ierr); 11724698041cSHong Zhang ierr = MatSetSizes(submats[0],0,0,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 11734698041cSHong Zhang ierr = MatSetType(submats[0],MATDUMMY);CHKERRQ(ierr); 11744698041cSHong Zhang 11755c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 11764698041cSHong Zhang ierr = PetscNewLog(submats[0],&smat_i);CHKERRQ(ierr); 11774698041cSHong Zhang submats[0]->data = (void*)smat_i; 11784698041cSHong Zhang 11794698041cSHong Zhang smat_i->destroy = submats[0]->ops->destroy; 1180f68bb481SHong Zhang submats[0]->ops->destroy = MatDestroySubMatrix_Dummy; 11814698041cSHong Zhang submats[0]->factortype = C->factortype; 11824698041cSHong Zhang 1183006cef31SHong Zhang smat_i->id = 0; 11844698041cSHong Zhang smat_i->nrqs = nrqs; 11854698041cSHong Zhang smat_i->nrqr = nrqr; 11864698041cSHong Zhang smat_i->rbuf1 = rbuf1; 11874698041cSHong Zhang smat_i->rbuf2 = rbuf2; 11884698041cSHong Zhang smat_i->rbuf3 = rbuf3; 11894698041cSHong Zhang smat_i->sbuf2 = sbuf2; 11904698041cSHong Zhang smat_i->req_source2 = req_source2; 11914698041cSHong Zhang 11924698041cSHong Zhang smat_i->sbuf1 = sbuf1; 11934698041cSHong Zhang smat_i->ptr = ptr; 11944698041cSHong Zhang smat_i->tmp = tmp; 11954698041cSHong Zhang smat_i->ctr = ctr; 11964698041cSHong Zhang 11974698041cSHong Zhang smat_i->pa = pa; 11984698041cSHong Zhang smat_i->req_size = req_size; 11994698041cSHong Zhang smat_i->req_source1 = req_source1; 12004698041cSHong Zhang 1201c67fe082SHong Zhang smat_i->allcolumns = PETSC_FALSE; 12024698041cSHong Zhang smat_i->singleis = PETSC_FALSE; 12034698041cSHong Zhang smat_i->row2proc = NULL; 12044698041cSHong Zhang smat_i->rmap = NULL; 12054698041cSHong Zhang smat_i->cmap = NULL; 12064698041cSHong Zhang } 12074698041cSHong Zhang 120817df9f7cSHong Zhang if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);} 120917df9f7cSHong Zhang ierr = PetscFree(lens);CHKERRQ(ierr); 1210175c2bc5SJunchao Zhang if (sbuf_aj) { 121117df9f7cSHong Zhang ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 121217df9f7cSHong Zhang ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1213175c2bc5SJunchao Zhang } 121417df9f7cSHong Zhang 121517df9f7cSHong Zhang } /* endof scall == MAT_INITIAL_MATRIX */ 121617df9f7cSHong Zhang 121717df9f7cSHong Zhang /* Post recv matrix values */ 1218bbc89d27SHong Zhang if (!ijonly) { 121917df9f7cSHong Zhang ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr); 1220175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&rbuf4);CHKERRQ(ierr); 1221175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqs,&r_waits4);CHKERRQ(ierr); 122217df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 1223a42ab0d6SHong Zhang ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr); 1224ffc4695bSBarry Smith ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRMPI(ierr); 122517df9f7cSHong Zhang } 122617df9f7cSHong Zhang 122717df9f7cSHong Zhang /* Allocate sending buffers for a->a, and send them off */ 1228175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&sbuf_aa);CHKERRQ(ierr); 122917df9f7cSHong Zhang for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 12309e686edcSHong Zhang 1231175c2bc5SJunchao Zhang if (nrqr) {ierr = PetscMalloc1(j*bs2,&sbuf_aa[0]);CHKERRQ(ierr);} 1232a42ab0d6SHong Zhang for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 123317df9f7cSHong Zhang 1234175c2bc5SJunchao Zhang ierr = PetscMalloc1(nrqr,&s_waits4);CHKERRQ(ierr); 123517df9f7cSHong Zhang 123617df9f7cSHong Zhang for (i=0; i<nrqr; i++) { 123717df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 123817df9f7cSHong Zhang sbuf_aa_i = sbuf_aa[i]; 123917df9f7cSHong Zhang ct1 = 2*rbuf1_i[0]+1; 124017df9f7cSHong Zhang ct2 = 0; 124117df9f7cSHong Zhang for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 124217df9f7cSHong Zhang kmax = rbuf1_i[2*j]; 124317df9f7cSHong Zhang for (k=0; k<kmax; k++,ct1++) { 124417df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 1245a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1246a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 124717df9f7cSHong Zhang ncols = nzA + nzB; 124817df9f7cSHong Zhang cworkB = b_j + b_i[row]; 1249a42ab0d6SHong Zhang vworkA = a_a + a_i[row]*bs2; 1250a42ab0d6SHong Zhang vworkB = b_a + b_i[row]*bs2; 125117df9f7cSHong Zhang 125217df9f7cSHong Zhang /* load the column values for this row into vals*/ 1253a42ab0d6SHong Zhang vals = sbuf_aa_i+ct2*bs2; 1254a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1255a42ab0d6SHong Zhang if ((bmap[cworkB[l]]) < cstart) { 1256580bdb30SBarry Smith ierr = PetscArraycpy(vals+l*bs2,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1257a42ab0d6SHong Zhang } else break; 1258a42ab0d6SHong Zhang } 1259a42ab0d6SHong Zhang imark = l; 1260a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1261580bdb30SBarry Smith ierr = PetscArraycpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2);CHKERRQ(ierr); 1262a42ab0d6SHong Zhang } 1263a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1264580bdb30SBarry Smith ierr = PetscArraycpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1265a42ab0d6SHong Zhang } 126617df9f7cSHong Zhang 126717df9f7cSHong Zhang ct2 += ncols; 126817df9f7cSHong Zhang } 126917df9f7cSHong Zhang } 1270ffc4695bSBarry Smith ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRMPI(ierr); 127117df9f7cSHong Zhang } 1272bbc89d27SHong Zhang } 127317df9f7cSHong Zhang 127417df9f7cSHong Zhang /* Assemble the matrices */ 127517df9f7cSHong Zhang /* First assemble the local rows */ 127617df9f7cSHong Zhang for (i=0; i<ismax; i++) { 127717df9f7cSHong Zhang row2proc_i = row2proc[i]; 127817df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 127917df9f7cSHong Zhang imat_ilen = subc->ilen; 128017df9f7cSHong Zhang imat_j = subc->j; 128117df9f7cSHong Zhang imat_i = subc->i; 128217df9f7cSHong Zhang imat_a = subc->a; 128317df9f7cSHong Zhang 128417df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 128517df9f7cSHong Zhang rmap_i = rmap[i]; 128617df9f7cSHong Zhang irow_i = irow[i]; 128717df9f7cSHong Zhang jmax = nrow[i]; 128817df9f7cSHong Zhang for (j=0; j<jmax; j++) { 1289ec3c739cSHong Zhang if (allrows[i]) row = j; 1290ec3c739cSHong Zhang else row = irow_i[j]; 129117df9f7cSHong Zhang proc = row2proc_i[j]; 1292a42ab0d6SHong Zhang 129317df9f7cSHong Zhang if (proc == rank) { 1294a42ab0d6SHong Zhang 1295a42ab0d6SHong Zhang row = row - rstart; 1296a42ab0d6SHong Zhang nzA = a_i[row+1] - a_i[row]; 1297a42ab0d6SHong Zhang nzB = b_i[row+1] - b_i[row]; 1298a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1299a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1300a42ab0d6SHong Zhang if (!ijonly) { 1301a42ab0d6SHong Zhang vworkA = a_a + a_i[row]*bs2; 1302a42ab0d6SHong Zhang vworkB = b_a + b_i[row]*bs2; 1303a42ab0d6SHong Zhang } 1304059f6b74SHong Zhang 1305059f6b74SHong Zhang if (allrows[i]) { 1306059f6b74SHong Zhang row = row+rstart; 1307059f6b74SHong Zhang } else { 1308a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1309a42ab0d6SHong Zhang ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1310a42ab0d6SHong Zhang row--; 1311121971b2SHong Zhang 1312a42ab0d6SHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1313a42ab0d6SHong Zhang #else 1314a42ab0d6SHong Zhang row = rmap_i[row + rstart]; 1315a42ab0d6SHong Zhang #endif 1316059f6b74SHong Zhang } 1317a42ab0d6SHong Zhang mat_i = imat_i[row]; 1318a42ab0d6SHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 1319a42ab0d6SHong Zhang mat_j = imat_j + mat_i; 132080d31651SHong Zhang ilen = imat_ilen[row]; 1321a42ab0d6SHong Zhang 1322a42ab0d6SHong Zhang /* load the column indices for this row into cols*/ 1323a42ab0d6SHong Zhang if (!allcolumns[i]) { 1324a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1325a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1326a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1327a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1328a42ab0d6SHong Zhang if (tcol) { 1329a42ab0d6SHong Zhang #else 1330a42ab0d6SHong Zhang if ((tcol = cmap_i[ctmp])) { 1331a42ab0d6SHong Zhang #endif 1332a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1333580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1334a42ab0d6SHong Zhang mat_a += bs2; 133580d31651SHong Zhang ilen++; 1336a42ab0d6SHong Zhang } 1337a42ab0d6SHong Zhang } else break; 1338a42ab0d6SHong Zhang } 1339a42ab0d6SHong Zhang imark = l; 1340a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1341a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1342a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1343a42ab0d6SHong Zhang if (tcol) { 1344a42ab0d6SHong Zhang #else 1345a42ab0d6SHong Zhang if ((tcol = cmap_i[cstart + cworkA[l]])) { 1346a42ab0d6SHong Zhang #endif 1347a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1348a42ab0d6SHong Zhang if (!ijonly) { 1349580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkA+l*bs2,bs2);CHKERRQ(ierr); 1350a42ab0d6SHong Zhang mat_a += bs2; 1351a42ab0d6SHong Zhang } 135280d31651SHong Zhang ilen++; 1353a42ab0d6SHong Zhang } 1354a42ab0d6SHong Zhang } 1355a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1356a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1357a42ab0d6SHong Zhang ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1358a42ab0d6SHong Zhang if (tcol) { 1359a42ab0d6SHong Zhang #else 1360a42ab0d6SHong Zhang if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1361a42ab0d6SHong Zhang #endif 1362a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1363a42ab0d6SHong Zhang if (!ijonly) { 1364580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1365a42ab0d6SHong Zhang mat_a += bs2; 1366a42ab0d6SHong Zhang } 136780d31651SHong Zhang ilen++; 1368a42ab0d6SHong Zhang } 1369a42ab0d6SHong Zhang } 1370a42ab0d6SHong Zhang } else { /* allcolumns */ 1371a42ab0d6SHong Zhang for (l=0; l<nzB; l++) { 1372a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1373a42ab0d6SHong Zhang *mat_j++ = ctmp; 1374580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1375a42ab0d6SHong Zhang mat_a += bs2; 137680d31651SHong Zhang ilen++; 1377a42ab0d6SHong Zhang } else break; 1378a42ab0d6SHong Zhang } 1379a42ab0d6SHong Zhang imark = l; 1380a42ab0d6SHong Zhang for (l=0; l<nzA; l++) { 1381a42ab0d6SHong Zhang *mat_j++ = cstart+cworkA[l]; 1382a42ab0d6SHong Zhang if (!ijonly) { 1383580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkA+l*bs2,bs2);CHKERRQ(ierr); 1384a42ab0d6SHong Zhang mat_a += bs2; 1385a42ab0d6SHong Zhang } 138680d31651SHong Zhang ilen++; 1387a42ab0d6SHong Zhang } 1388a42ab0d6SHong Zhang for (l=imark; l<nzB; l++) { 1389a42ab0d6SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1390a42ab0d6SHong Zhang if (!ijonly) { 1391580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,vworkB+l*bs2,bs2);CHKERRQ(ierr); 1392a42ab0d6SHong Zhang mat_a += bs2; 1393a42ab0d6SHong Zhang } 139480d31651SHong Zhang ilen++; 1395a42ab0d6SHong Zhang } 1396a42ab0d6SHong Zhang } 139780d31651SHong Zhang imat_ilen[row] = ilen; 139817df9f7cSHong Zhang } 139917df9f7cSHong Zhang } 140017df9f7cSHong Zhang } 140117df9f7cSHong Zhang 140217df9f7cSHong Zhang /* Now assemble the off proc rows */ 14035dd0c08cSHong Zhang if (!ijonly) { 1404175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqs,r_waits4,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 14055dd0c08cSHong Zhang } 140617df9f7cSHong Zhang for (tmp2=0; tmp2<nrqs; tmp2++) { 140717df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 140817df9f7cSHong Zhang jmax = sbuf1_i[0]; 140917df9f7cSHong Zhang ct1 = 2*jmax + 1; 141017df9f7cSHong Zhang ct2 = 0; 141117df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 141217df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 14135dd0c08cSHong Zhang if (!ijonly) rbuf4_i = rbuf4[tmp2]; 141417df9f7cSHong Zhang for (j=1; j<=jmax; j++) { 141517df9f7cSHong Zhang is_no = sbuf1_i[2*j-1]; 141617df9f7cSHong Zhang rmap_i = rmap[is_no]; 141717df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 141817df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[is_no]->data; 141917df9f7cSHong Zhang imat_ilen = subc->ilen; 142017df9f7cSHong Zhang imat_j = subc->j; 142117df9f7cSHong Zhang imat_i = subc->i; 14225dd0c08cSHong Zhang if (!ijonly) imat_a = subc->a; 142317df9f7cSHong Zhang max1 = sbuf1_i[2*j]; 14245dd0c08cSHong Zhang for (k=0; k<max1; k++,ct1++) { /* for each recved block row */ 142517df9f7cSHong Zhang row = sbuf1_i[ct1]; 1426059f6b74SHong Zhang 1427059f6b74SHong Zhang if (allrows[is_no]) { 1428059f6b74SHong Zhang row = sbuf1_i[ct1]; 1429059f6b74SHong Zhang } else { 143017df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 143117df9f7cSHong Zhang ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 143217df9f7cSHong Zhang row--; 14335dd0c08cSHong Zhang if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 143417df9f7cSHong Zhang #else 143517df9f7cSHong Zhang row = rmap_i[row]; 143617df9f7cSHong Zhang #endif 1437059f6b74SHong Zhang } 143817df9f7cSHong Zhang ilen = imat_ilen[row]; 143917df9f7cSHong Zhang mat_i = imat_i[row]; 14405dd0c08cSHong Zhang if (!ijonly) mat_a = imat_a + mat_i*bs2; 144117df9f7cSHong Zhang mat_j = imat_j + mat_i; 144217df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 144317df9f7cSHong Zhang if (!allcolumns[is_no]) { 144417df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 144517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 144617df9f7cSHong Zhang ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 144717df9f7cSHong Zhang #else 144817df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 144917df9f7cSHong Zhang #endif 145017df9f7cSHong Zhang if (tcol) { 145117df9f7cSHong Zhang *mat_j++ = tcol - 1; 14525dd0c08cSHong Zhang if (!ijonly) { 1453580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2);CHKERRQ(ierr); 14545dd0c08cSHong Zhang mat_a += bs2; 14555dd0c08cSHong Zhang } 145617df9f7cSHong Zhang ilen++; 145717df9f7cSHong Zhang } 145817df9f7cSHong Zhang } 145917df9f7cSHong Zhang } else { /* allcolumns */ 146017df9f7cSHong Zhang for (l=0; l<max2; l++,ct2++) { 146117df9f7cSHong Zhang *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 14625dd0c08cSHong Zhang if (!ijonly) { 1463580bdb30SBarry Smith ierr = PetscArraycpy(mat_a,rbuf4_i+ct2*bs2,bs2);CHKERRQ(ierr); 14645dd0c08cSHong Zhang mat_a += bs2; 14655dd0c08cSHong Zhang } 146617df9f7cSHong Zhang ilen++; 146717df9f7cSHong Zhang } 146817df9f7cSHong Zhang } 146917df9f7cSHong Zhang imat_ilen[row] = ilen; 147017df9f7cSHong Zhang } 147117df9f7cSHong Zhang } 147217df9f7cSHong Zhang } 147317df9f7cSHong Zhang 147480d31651SHong Zhang if (!iscsorted) { /* sort column indices of the rows */ 147580d31651SHong Zhang MatScalar *work; 147680d31651SHong Zhang 147780d31651SHong Zhang ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 147817df9f7cSHong Zhang for (i=0; i<ismax; i++) { 147917df9f7cSHong Zhang subc = (Mat_SeqBAIJ*)submats[i]->data; 148080d31651SHong Zhang imat_ilen = subc->ilen; 148117df9f7cSHong Zhang imat_j = subc->j; 148217df9f7cSHong Zhang imat_i = subc->i; 14834b6ceb0dSHong Zhang if (!ijonly) imat_a = subc->a; 148417df9f7cSHong Zhang if (allcolumns[i]) continue; 14854b6ceb0dSHong Zhang 148617df9f7cSHong Zhang jmax = nrow[i]; 148717df9f7cSHong Zhang for (j=0; j<jmax; j++) { 148880d31651SHong Zhang mat_i = imat_i[j]; 148980d31651SHong Zhang mat_j = imat_j + mat_i; 14904b6ceb0dSHong Zhang ilen = imat_ilen[j]; 149180d31651SHong Zhang if (ijonly) { 149280d31651SHong Zhang ierr = PetscSortInt(ilen,mat_j);CHKERRQ(ierr); 149380d31651SHong Zhang } else { 14944b6ceb0dSHong Zhang mat_a = imat_a + mat_i*bs2; 149580d31651SHong Zhang ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr); 149617df9f7cSHong Zhang } 149717df9f7cSHong Zhang } 149817df9f7cSHong Zhang } 149980d31651SHong Zhang ierr = PetscFree(work);CHKERRQ(ierr); 150080d31651SHong Zhang } 150117df9f7cSHong Zhang 1502bbc89d27SHong Zhang if (!ijonly) { 150317df9f7cSHong Zhang ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1504175c2bc5SJunchao Zhang ierr = MPI_Waitall(nrqr,s_waits4,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 150517df9f7cSHong Zhang ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1506bbc89d27SHong Zhang } 150717df9f7cSHong Zhang 150817df9f7cSHong Zhang /* Restore the indices */ 150917df9f7cSHong Zhang for (i=0; i<ismax; i++) { 1510ec3c739cSHong Zhang if (!allrows[i]) { 151117df9f7cSHong Zhang ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1512ec3c739cSHong Zhang } 151317df9f7cSHong Zhang if (!allcolumns[i]) { 151417df9f7cSHong Zhang ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 151517df9f7cSHong Zhang } 151617df9f7cSHong Zhang } 151717df9f7cSHong Zhang 151817df9f7cSHong Zhang for (i=0; i<ismax; i++) { 151917df9f7cSHong Zhang ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152017df9f7cSHong Zhang ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152117df9f7cSHong Zhang } 152217df9f7cSHong Zhang 15232cd485c2SBarry Smith ierr = PetscFree5(*(PetscInt***)&irow,*(PetscInt***)&icol,nrow,ncol,issorted);CHKERRQ(ierr); 15244698041cSHong Zhang ierr = PetscFree5(row2proc,cmap,rmap,allcolumns,allrows);CHKERRQ(ierr); 152517df9f7cSHong Zhang 1526bbc89d27SHong Zhang if (!ijonly) { 1527175c2bc5SJunchao Zhang if (sbuf_aa) { 152817df9f7cSHong Zhang ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 152917df9f7cSHong Zhang ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1530175c2bc5SJunchao Zhang } 153117df9f7cSHong Zhang 153217df9f7cSHong Zhang for (i=0; i<nrqs; ++i) { 153317df9f7cSHong Zhang ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 153417df9f7cSHong Zhang } 153517df9f7cSHong Zhang ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1536bbc89d27SHong Zhang } 15377a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 15383a40ed3dSBarry Smith PetscFunctionReturn(0); 1539d5b485abSSatish Balay } 1540