18965ea79SLois Curfman McInnes #ifndef lint 2*3501a2bdSLois Curfman McInnes static char vcid[] = "$Id: mpidense.c,v 1.4 1995/10/22 22:17:20 curfman Exp curfman $"; 38965ea79SLois Curfman McInnes #endif 48965ea79SLois Curfman McInnes 58965ea79SLois Curfman McInnes #include "mpidense.h" 68965ea79SLois Curfman McInnes #include "vec/vecimpl.h" 78965ea79SLois Curfman McInnes #include "inline/spops.h" 88965ea79SLois Curfman McInnes 939ddd567SLois Curfman McInnes static int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n, 108965ea79SLois Curfman McInnes int *idxn,Scalar *v,InsertMode addv) 118965ea79SLois Curfman McInnes { 1239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1339ddd567SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 148965ea79SLois Curfman McInnes 1539ddd567SLois Curfman McInnes if (mdn->insertmode != NOT_SET_VALUES && mdn->insertmode != addv) { 1639ddd567SLois Curfman McInnes SETERRQ(1,"MatSetValues_MPIDense:Cannot mix inserts and adds"); 178965ea79SLois Curfman McInnes } 1839ddd567SLois Curfman McInnes mdn->insertmode = addv; 198965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 2039ddd567SLois Curfman McInnes if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative row"); 2139ddd567SLois Curfman McInnes if (idxm[i] >= mdn->M) SETERRQ(1,"MatSetValues_MPIDense:Row too large"); 228965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 238965ea79SLois Curfman McInnes row = idxm[i] - rstart; 248965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 2539ddd567SLois Curfman McInnes if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative column"); 2639ddd567SLois Curfman McInnes if (idxn[j] >= mdn->N) 2739ddd567SLois Curfman McInnes SETERRQ(1,"MatSetValues_MPIDense:Column too large"); 2839ddd567SLois Curfman McInnes ierr = MatSetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j,addv); 2939ddd567SLois Curfman McInnes CHKERRQ(ierr); 308965ea79SLois Curfman McInnes } 318965ea79SLois Curfman McInnes } 328965ea79SLois Curfman McInnes else { 3339ddd567SLois Curfman McInnes ierr = StashValues_Private(&mdn->stash,idxm[i],n,idxn,v+i*n,addv); 3439ddd567SLois Curfman McInnes CHKERRQ(ierr); 358965ea79SLois Curfman McInnes } 368965ea79SLois Curfman McInnes } 378965ea79SLois Curfman McInnes return 0; 388965ea79SLois Curfman McInnes } 398965ea79SLois Curfman McInnes 4039ddd567SLois Curfman McInnes static int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 418965ea79SLois Curfman McInnes { 4239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 438965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 4439ddd567SLois Curfman McInnes int size = mdn->size, *owners = mdn->rowners, rank = mdn->rank; 458965ea79SLois Curfman McInnes int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 4639ddd567SLois Curfman McInnes int tag = mat->tag, *owner,*starts,count,ierr; 478965ea79SLois Curfman McInnes InsertMode addv; 4839ddd567SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 498965ea79SLois Curfman McInnes Scalar *rvalues,*svalues; 508965ea79SLois Curfman McInnes 518965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 5239ddd567SLois Curfman McInnes MPI_Allreduce((void *) &mdn->insertmode,(void *) &addv,1,MPI_INT, 5339ddd567SLois Curfman McInnes MPI_BOR,comm); 5439ddd567SLois Curfman McInnes if (addv == (ADD_VALUES|INSERT_VALUES)) { SETERRQ(1, 5539ddd567SLois Curfman McInnes "MatAssemblyBegin_MPIDense:Cannot mix adds/inserts on different procs"); 568965ea79SLois Curfman McInnes } 5739ddd567SLois Curfman McInnes mdn->insertmode = addv; /* in case this processor had no cache */ 588965ea79SLois Curfman McInnes 598965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 608965ea79SLois Curfman McInnes nprocs = (int *) PETSCMALLOC( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 618965ea79SLois Curfman McInnes PetscZero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 6239ddd567SLois Curfman McInnes owner = (int *) PETSCMALLOC( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 6339ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 6439ddd567SLois Curfman McInnes idx = mdn->stash.idx[i]; 658965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 668965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 678965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; break; 688965ea79SLois Curfman McInnes } 698965ea79SLois Curfman McInnes } 708965ea79SLois Curfman McInnes } 718965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 728965ea79SLois Curfman McInnes 738965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 748965ea79SLois Curfman McInnes work = (int *) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(work); 7539ddd567SLois Curfman McInnes MPI_Allreduce((void *) procs,(void *) work,size,MPI_INT,MPI_SUM,comm); 768965ea79SLois Curfman McInnes nreceives = work[rank]; 7739ddd567SLois Curfman McInnes MPI_Allreduce((void *) nprocs,(void *) work,size,MPI_INT,MPI_MAX,comm); 788965ea79SLois Curfman McInnes nmax = work[rank]; 798965ea79SLois Curfman McInnes PETSCFREE(work); 808965ea79SLois Curfman McInnes 818965ea79SLois Curfman McInnes /* post receives: 828965ea79SLois Curfman McInnes 1) each message will consist of ordered pairs 838965ea79SLois Curfman McInnes (global index,value) we store the global index as a double 848965ea79SLois Curfman McInnes to simplify the message passing. 858965ea79SLois Curfman McInnes 2) since we don't know how long each individual message is we 868965ea79SLois Curfman McInnes allocate the largest needed buffer for each receive. Potentially 878965ea79SLois Curfman McInnes this is a lot of wasted space. 888965ea79SLois Curfman McInnes 898965ea79SLois Curfman McInnes This could be done better. 908965ea79SLois Curfman McInnes */ 918965ea79SLois Curfman McInnes rvalues = (Scalar *) PETSCMALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 928965ea79SLois Curfman McInnes CHKPTRQ(rvalues); 938965ea79SLois Curfman McInnes recv_waits = (MPI_Request *) PETSCMALLOC((nreceives+1)*sizeof(MPI_Request)); 948965ea79SLois Curfman McInnes CHKPTRQ(recv_waits); 958965ea79SLois Curfman McInnes for ( i=0; i<nreceives; i++ ) { 9639ddd567SLois Curfman McInnes MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 978965ea79SLois Curfman McInnes comm,recv_waits+i); 988965ea79SLois Curfman McInnes } 998965ea79SLois Curfman McInnes 1008965ea79SLois Curfman McInnes /* do sends: 1018965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 1028965ea79SLois Curfman McInnes the ith processor 1038965ea79SLois Curfman McInnes */ 10439ddd567SLois Curfman McInnes svalues = (Scalar *) PETSCMALLOC( 3*(mdn->stash.n+1)*sizeof(Scalar) ); 10539ddd567SLois Curfman McInnes CHKPTRQ(svalues); 1068965ea79SLois Curfman McInnes send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request)); 1078965ea79SLois Curfman McInnes CHKPTRQ(send_waits); 1088965ea79SLois Curfman McInnes starts = (int *) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(starts); 1098965ea79SLois Curfman McInnes starts[0] = 0; 1108965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 11139ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 11239ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]] = (Scalar) mdn->stash.idx[i]; 11339ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]+1] = (Scalar) mdn->stash.idy[i]; 11439ddd567SLois Curfman McInnes svalues[3*(starts[owner[i]]++)+2] = mdn->stash.array[i]; 1158965ea79SLois Curfman McInnes } 1168965ea79SLois Curfman McInnes PETSCFREE(owner); 1178965ea79SLois Curfman McInnes starts[0] = 0; 1188965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1198965ea79SLois Curfman McInnes count = 0; 1208965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 1218965ea79SLois Curfman McInnes if (procs[i]) { 12239ddd567SLois Curfman McInnes MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPIU_SCALAR,i,tag, 1238965ea79SLois Curfman McInnes comm,send_waits+count++); 1248965ea79SLois Curfman McInnes } 1258965ea79SLois Curfman McInnes } 1268965ea79SLois Curfman McInnes PETSCFREE(starts); PETSCFREE(nprocs); 1278965ea79SLois Curfman McInnes 1288965ea79SLois Curfman McInnes /* Free cache space */ 12939ddd567SLois Curfman McInnes ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr); 1308965ea79SLois Curfman McInnes 13139ddd567SLois Curfman McInnes mdn->svalues = svalues; mdn->rvalues = rvalues; 13239ddd567SLois Curfman McInnes mdn->nsends = nsends; mdn->nrecvs = nreceives; 13339ddd567SLois Curfman McInnes mdn->send_waits = send_waits; mdn->recv_waits = recv_waits; 13439ddd567SLois Curfman McInnes mdn->rmax = nmax; 1358965ea79SLois Curfman McInnes 1368965ea79SLois Curfman McInnes return 0; 1378965ea79SLois Curfman McInnes } 13839ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat); 1398965ea79SLois Curfman McInnes 14039ddd567SLois Curfman McInnes static int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1418965ea79SLois Curfman McInnes { 14239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1438965ea79SLois Curfman McInnes MPI_Status *send_status,recv_status; 14439ddd567SLois Curfman McInnes int imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col; 1458965ea79SLois Curfman McInnes Scalar *values,val; 14639ddd567SLois Curfman McInnes InsertMode addv = mdn->insertmode; 1478965ea79SLois Curfman McInnes 1488965ea79SLois Curfman McInnes /* wait on receives */ 1498965ea79SLois Curfman McInnes while (count) { 15039ddd567SLois Curfman McInnes MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status); 1518965ea79SLois Curfman McInnes /* unpack receives into our local space */ 15239ddd567SLois Curfman McInnes values = mdn->rvalues + 3*imdex*mdn->rmax; 1538965ea79SLois Curfman McInnes MPI_Get_count(&recv_status,MPIU_SCALAR,&n); 1548965ea79SLois Curfman McInnes n = n/3; 1558965ea79SLois Curfman McInnes for ( i=0; i<n; i++ ) { 15639ddd567SLois Curfman McInnes row = (int) PETSCREAL(values[3*i]) - mdn->rstart; 1578965ea79SLois Curfman McInnes col = (int) PETSCREAL(values[3*i+1]); 1588965ea79SLois Curfman McInnes val = values[3*i+2]; 15939ddd567SLois Curfman McInnes if (col >= 0 && col < mdn->N) { 16039ddd567SLois Curfman McInnes MatSetValues(mdn->A,1,&row,1,&col,&val,addv); 1618965ea79SLois Curfman McInnes } 16239ddd567SLois Curfman McInnes else {SETERRQ(1,"MatAssemblyEnd_MPIDense:Invalid column");} 1638965ea79SLois Curfman McInnes } 1648965ea79SLois Curfman McInnes count--; 1658965ea79SLois Curfman McInnes } 16639ddd567SLois Curfman McInnes PETSCFREE(mdn->recv_waits); PETSCFREE(mdn->rvalues); 1678965ea79SLois Curfman McInnes 1688965ea79SLois Curfman McInnes /* wait on sends */ 16939ddd567SLois Curfman McInnes if (mdn->nsends) { 17039ddd567SLois Curfman McInnes send_status = (MPI_Status *) PETSCMALLOC( mdn->nsends*sizeof(MPI_Status) ); 1718965ea79SLois Curfman McInnes CHKPTRQ(send_status); 17239ddd567SLois Curfman McInnes MPI_Waitall(mdn->nsends,mdn->send_waits,send_status); 1738965ea79SLois Curfman McInnes PETSCFREE(send_status); 1748965ea79SLois Curfman McInnes } 17539ddd567SLois Curfman McInnes PETSCFREE(mdn->send_waits); PETSCFREE(mdn->svalues); 1768965ea79SLois Curfman McInnes 17739ddd567SLois Curfman McInnes mdn->insertmode = NOT_SET_VALUES; 17839ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr); 17939ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr); 1808965ea79SLois Curfman McInnes 18139ddd567SLois Curfman McInnes if (!mdn->assembled && mode == FINAL_ASSEMBLY) { 18239ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr); 1838965ea79SLois Curfman McInnes } 18439ddd567SLois Curfman McInnes mdn->assembled = 1; 1858965ea79SLois Curfman McInnes return 0; 1868965ea79SLois Curfman McInnes } 1878965ea79SLois Curfman McInnes 18839ddd567SLois Curfman McInnes static int MatZeroEntries_MPIDense(Mat A) 1898965ea79SLois Curfman McInnes { 19039ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 19139ddd567SLois Curfman McInnes return MatZeroEntries(l->A); 1928965ea79SLois Curfman McInnes } 1938965ea79SLois Curfman McInnes 1948965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 1958965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 1968965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 197*3501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 1988965ea79SLois Curfman McInnes routine. 1998965ea79SLois Curfman McInnes */ 20039ddd567SLois Curfman McInnes static int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 2018965ea79SLois Curfman McInnes { 20239ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2038965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 2048965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 2058965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2068965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2078965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2088965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2098965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2108965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2118965ea79SLois Curfman McInnes IS istmp; 2128965ea79SLois Curfman McInnes 21339ddd567SLois Curfman McInnes if (!l->assembled) SETERRQ(1,"MatZeroRows_MPIDense:Must assemble matrix"); 2148965ea79SLois Curfman McInnes ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 2158965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 2168965ea79SLois Curfman McInnes 2178965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2188965ea79SLois Curfman McInnes nprocs = (int *) PETSCMALLOC( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 2198965ea79SLois Curfman McInnes PetscZero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 2208965ea79SLois Curfman McInnes owner = (int *) PETSCMALLOC((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 2218965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2228965ea79SLois Curfman McInnes idx = rows[i]; 2238965ea79SLois Curfman McInnes found = 0; 2248965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2258965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2268965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2278965ea79SLois Curfman McInnes } 2288965ea79SLois Curfman McInnes } 22939ddd567SLois Curfman McInnes if (!found) SETERRQ(1,"MatZeroRows_MPIDense:Index out of range"); 2308965ea79SLois Curfman McInnes } 2318965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 2328965ea79SLois Curfman McInnes 2338965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 2348965ea79SLois Curfman McInnes work = (int *) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(work); 2358965ea79SLois Curfman McInnes MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm); 2368965ea79SLois Curfman McInnes nrecvs = work[rank]; 2378965ea79SLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 2388965ea79SLois Curfman McInnes nmax = work[rank]; 2398965ea79SLois Curfman McInnes PETSCFREE(work); 2408965ea79SLois Curfman McInnes 2418965ea79SLois Curfman McInnes /* post receives: */ 2428965ea79SLois Curfman McInnes rvalues = (int *) PETSCMALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 2438965ea79SLois Curfman McInnes CHKPTRQ(rvalues); 2448965ea79SLois Curfman McInnes recv_waits = (MPI_Request *) PETSCMALLOC((nrecvs+1)*sizeof(MPI_Request)); 2458965ea79SLois Curfman McInnes CHKPTRQ(recv_waits); 2468965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2478965ea79SLois Curfman McInnes MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i); 2488965ea79SLois Curfman McInnes } 2498965ea79SLois Curfman McInnes 2508965ea79SLois Curfman McInnes /* do sends: 2518965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 2528965ea79SLois Curfman McInnes the ith processor 2538965ea79SLois Curfman McInnes */ 2548965ea79SLois Curfman McInnes svalues = (int *) PETSCMALLOC( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 2558965ea79SLois Curfman McInnes send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request)); 2568965ea79SLois Curfman McInnes CHKPTRQ(send_waits); 2578965ea79SLois Curfman McInnes starts = (int *) PETSCMALLOC( (size+1)*sizeof(int) ); CHKPTRQ(starts); 2588965ea79SLois Curfman McInnes starts[0] = 0; 2598965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2608965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2618965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 2628965ea79SLois Curfman McInnes } 2638965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 2648965ea79SLois Curfman McInnes 2658965ea79SLois Curfman McInnes starts[0] = 0; 2668965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2678965ea79SLois Curfman McInnes count = 0; 2688965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 2698965ea79SLois Curfman McInnes if (procs[i]) { 2708965ea79SLois Curfman McInnes MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++); 2718965ea79SLois Curfman McInnes } 2728965ea79SLois Curfman McInnes } 2738965ea79SLois Curfman McInnes PETSCFREE(starts); 2748965ea79SLois Curfman McInnes 2758965ea79SLois Curfman McInnes base = owners[rank]; 2768965ea79SLois Curfman McInnes 2778965ea79SLois Curfman McInnes /* wait on receives */ 2788965ea79SLois Curfman McInnes lens = (int *) PETSCMALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 2798965ea79SLois Curfman McInnes source = lens + nrecvs; 2808965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 2818965ea79SLois Curfman McInnes while (count) { 2828965ea79SLois Curfman McInnes MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 2838965ea79SLois Curfman McInnes /* unpack receives into our local space */ 2848965ea79SLois Curfman McInnes MPI_Get_count(&recv_status,MPI_INT,&n); 2858965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 2868965ea79SLois Curfman McInnes lens[imdex] = n; 2878965ea79SLois Curfman McInnes slen += n; 2888965ea79SLois Curfman McInnes count--; 2898965ea79SLois Curfman McInnes } 2908965ea79SLois Curfman McInnes PETSCFREE(recv_waits); 2918965ea79SLois Curfman McInnes 2928965ea79SLois Curfman McInnes /* move the data into the send scatter */ 2938965ea79SLois Curfman McInnes lrows = (int *) PETSCMALLOC( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 2948965ea79SLois Curfman McInnes count = 0; 2958965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2968965ea79SLois Curfman McInnes values = rvalues + i*nmax; 2978965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 2988965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 2998965ea79SLois Curfman McInnes } 3008965ea79SLois Curfman McInnes } 3018965ea79SLois Curfman McInnes PETSCFREE(rvalues); PETSCFREE(lens); 3028965ea79SLois Curfman McInnes PETSCFREE(owner); PETSCFREE(nprocs); 3038965ea79SLois Curfman McInnes 3048965ea79SLois Curfman McInnes /* actually zap the local rows */ 3058965ea79SLois Curfman McInnes ierr = ISCreateSeq(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 3068965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 3078965ea79SLois Curfman McInnes PETSCFREE(lrows); 3088965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 3098965ea79SLois Curfman McInnes ierr = ISDestroy(istmp); CHKERRQ(ierr); 3108965ea79SLois Curfman McInnes 3118965ea79SLois Curfman McInnes /* wait on sends */ 3128965ea79SLois Curfman McInnes if (nsends) { 3138965ea79SLois Curfman McInnes send_status = (MPI_Status *) PETSCMALLOC(nsends*sizeof(MPI_Status)); 3148965ea79SLois Curfman McInnes CHKPTRQ(send_status); 3158965ea79SLois Curfman McInnes MPI_Waitall(nsends,send_waits,send_status); 3168965ea79SLois Curfman McInnes PETSCFREE(send_status); 3178965ea79SLois Curfman McInnes } 3188965ea79SLois Curfman McInnes PETSCFREE(send_waits); PETSCFREE(svalues); 3198965ea79SLois Curfman McInnes 3208965ea79SLois Curfman McInnes return 0; 3218965ea79SLois Curfman McInnes } 3228965ea79SLois Curfman McInnes 32339ddd567SLois Curfman McInnes static int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3248965ea79SLois Curfman McInnes { 32539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3268965ea79SLois Curfman McInnes int ierr; 32739ddd567SLois Curfman McInnes if (!mdn->assembled) 32839ddd567SLois Curfman McInnes SETERRQ(1,"MatMult_MPIDense:Must assemble matrix first"); 32939ddd567SLois Curfman McInnes ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 33039ddd567SLois Curfman McInnes CHKERRQ(ierr); 33139ddd567SLois Curfman McInnes ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 33239ddd567SLois Curfman McInnes CHKERRQ(ierr); 33339ddd567SLois Curfman McInnes ierr = MatMult(mdn->A,mdn->lvec,yy); CHKERRQ(ierr); 3348965ea79SLois Curfman McInnes return 0; 3358965ea79SLois Curfman McInnes } 3368965ea79SLois Curfman McInnes 33739ddd567SLois Curfman McInnes static int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3388965ea79SLois Curfman McInnes { 33939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3408965ea79SLois Curfman McInnes int ierr; 34139ddd567SLois Curfman McInnes if (!mdn->assembled) 34239ddd567SLois Curfman McInnes SETERRQ(1,"MatMultAdd_MPIDense:Must assemble matrix first"); 343*3501a2bdSLois Curfman McInnes ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 34439ddd567SLois Curfman McInnes CHKERRQ(ierr); 345*3501a2bdSLois Curfman McInnes ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 34639ddd567SLois Curfman McInnes CHKERRQ(ierr); 34739ddd567SLois Curfman McInnes ierr = MatMultAdd(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr); 3488965ea79SLois Curfman McInnes return 0; 3498965ea79SLois Curfman McInnes } 3508965ea79SLois Curfman McInnes 351096963f5SLois Curfman McInnes static int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 352096963f5SLois Curfman McInnes { 353096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 354096963f5SLois Curfman McInnes int ierr; 355*3501a2bdSLois Curfman McInnes Scalar zero = 0.0; 356096963f5SLois Curfman McInnes 357096963f5SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatMulTrans_MPIDense:must assemble matrix"); 358*3501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy); CHKERRQ(ierr); 359096963f5SLois Curfman McInnes ierr = MatMultTrans(a->A,xx,a->lvec); CHKERRQ(ierr); 360096963f5SLois Curfman McInnes ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES, 361096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 362096963f5SLois Curfman McInnes ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES, 363096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 364096963f5SLois Curfman McInnes return 0; 365096963f5SLois Curfman McInnes } 366096963f5SLois Curfman McInnes 367096963f5SLois Curfman McInnes static int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 368096963f5SLois Curfman McInnes { 369096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 370096963f5SLois Curfman McInnes int ierr; 371096963f5SLois Curfman McInnes 372096963f5SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatMulTransAdd_MPIDense:must assemble matrix"); 373*3501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz); CHKERRQ(ierr); 374*3501a2bdSLois Curfman McInnes ierr = MatMultTrans(a->A,xx,a->lvec); CHKERRQ(ierr); 375096963f5SLois Curfman McInnes ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES, 376096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 377096963f5SLois Curfman McInnes ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES, 378096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 379096963f5SLois Curfman McInnes return 0; 380096963f5SLois Curfman McInnes } 381096963f5SLois Curfman McInnes 38239ddd567SLois Curfman McInnes static int MatGetDiagonal_MPIDense(Mat A,Vec v) 3838965ea79SLois Curfman McInnes { 38439ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 385096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 386096963f5SLois Curfman McInnes int ierr, i, n, m = a->m, radd; 387096963f5SLois Curfman McInnes Scalar *x; 38839ddd567SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatGetDiag_MPIDense:must assemble matrix"); 389096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x); CHKERRQ(ierr); 390096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n); CHKERRQ(ierr); 391096963f5SLois Curfman McInnes if (n != a->M) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 392096963f5SLois Curfman McInnes radd = a->rstart*m*m; 393096963f5SLois Curfman McInnes for ( i=0; i<m; i++ ) { 394096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 395096963f5SLois Curfman McInnes } 396096963f5SLois Curfman McInnes return 0; 3978965ea79SLois Curfman McInnes } 3988965ea79SLois Curfman McInnes 39939ddd567SLois Curfman McInnes static int MatDestroy_MPIDense(PetscObject obj) 4008965ea79SLois Curfman McInnes { 4018965ea79SLois Curfman McInnes Mat mat = (Mat) obj; 402*3501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4038965ea79SLois Curfman McInnes int ierr; 4048965ea79SLois Curfman McInnes #if defined(PETSC_LOG) 405*3501a2bdSLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4068965ea79SLois Curfman McInnes #endif 407*3501a2bdSLois Curfman McInnes PETSCFREE(mdn->rowners); 408*3501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A); CHKERRQ(ierr); 409*3501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 410*3501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 411*3501a2bdSLois Curfman McInnes PETSCFREE(mdn); 4128965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4138965ea79SLois Curfman McInnes PETSCHEADERDESTROY(mat); 4148965ea79SLois Curfman McInnes return 0; 4158965ea79SLois Curfman McInnes } 41639ddd567SLois Curfman McInnes 4178965ea79SLois Curfman McInnes #include "pinclude/pviewer.h" 4188965ea79SLois Curfman McInnes 41939ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4208965ea79SLois Curfman McInnes { 42139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4228965ea79SLois Curfman McInnes int ierr; 42339ddd567SLois Curfman McInnes if (mdn->size == 1) { 42439ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4258965ea79SLois Curfman McInnes } 42639ddd567SLois Curfman McInnes else SETERRQ(1,"MatView_MPIDense_Binary:Only uniprocessor output supported"); 4278965ea79SLois Curfman McInnes return 0; 4288965ea79SLois Curfman McInnes } 4298965ea79SLois Curfman McInnes 43039ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 4318965ea79SLois Curfman McInnes { 43239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 43339ddd567SLois Curfman McInnes int ierr, format; 4348965ea79SLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 4358965ea79SLois Curfman McInnes FILE *fd; 4368965ea79SLois Curfman McInnes 437*3501a2bdSLois Curfman McInnes ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 4388965ea79SLois Curfman McInnes if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 4398965ea79SLois Curfman McInnes ierr = ViewerFileGetFormat_Private(viewer,&format); 440*3501a2bdSLois Curfman McInnes if (format == FILE_FORMAT_INFO_DETAILED) { 441096963f5SLois Curfman McInnes int nz, nzalloc, mem, rank; 442096963f5SLois Curfman McInnes MPI_Comm_rank(mat->comm,&rank); 443096963f5SLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); 444096963f5SLois Curfman McInnes MPIU_Seq_begin(mat->comm,1); 445*3501a2bdSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n", 446096963f5SLois Curfman McInnes rank,mdn->m,nz,nzalloc,mem); 447096963f5SLois Curfman McInnes fflush(fd); 448096963f5SLois Curfman McInnes MPIU_Seq_end(mat->comm,1); 449*3501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr); 450*3501a2bdSLois Curfman McInnes return 0; 451*3501a2bdSLois Curfman McInnes } 452*3501a2bdSLois Curfman McInnes else if (format == FILE_FORMAT_INFO) { 453096963f5SLois Curfman McInnes return 0; 4548965ea79SLois Curfman McInnes } 4558965ea79SLois Curfman McInnes } 4568965ea79SLois Curfman McInnes if (vobj->type == ASCII_FILE_VIEWER) { 4578965ea79SLois Curfman McInnes MPIU_Seq_begin(mat->comm,1); 45839ddd567SLois Curfman McInnes fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n", 45939ddd567SLois Curfman McInnes mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n); 46039ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4618965ea79SLois Curfman McInnes fflush(fd); 4628965ea79SLois Curfman McInnes MPIU_Seq_end(mat->comm,1); 4638965ea79SLois Curfman McInnes } 4648965ea79SLois Curfman McInnes else { 46539ddd567SLois Curfman McInnes int size = mdn->size, rank = mdn->rank; 4668965ea79SLois Curfman McInnes if (size == 1) { 46739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4688965ea79SLois Curfman McInnes } 4698965ea79SLois Curfman McInnes else { 4708965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 4718965ea79SLois Curfman McInnes Mat A; 47239ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 47339ddd567SLois Curfman McInnes Scalar *vals; 47439ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 4758965ea79SLois Curfman McInnes 4768965ea79SLois Curfman McInnes if (!rank) { 47739ddd567SLois Curfman McInnes ierr = MatCreateMPIDense(mat->comm,M,M,N,N,&A); CHKERRQ(ierr); 4788965ea79SLois Curfman McInnes } 4798965ea79SLois Curfman McInnes else { 48039ddd567SLois Curfman McInnes ierr = MatCreateMPIDense(mat->comm,0,M,N,N,&A); CHKERRQ(ierr); 4818965ea79SLois Curfman McInnes } 4828965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 4838965ea79SLois Curfman McInnes 48439ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 48539ddd567SLois Curfman McInnes but it's quick for now */ 48639ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 4878965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 48839ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 48939ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 49039ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 49139ddd567SLois Curfman McInnes row++; 4928965ea79SLois Curfman McInnes } 4938965ea79SLois Curfman McInnes 4948965ea79SLois Curfman McInnes ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 4958965ea79SLois Curfman McInnes ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 4968965ea79SLois Curfman McInnes if (!rank) { 49739ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr); 4988965ea79SLois Curfman McInnes } 4998965ea79SLois Curfman McInnes ierr = MatDestroy(A); CHKERRQ(ierr); 5008965ea79SLois Curfman McInnes } 5018965ea79SLois Curfman McInnes } 5028965ea79SLois Curfman McInnes return 0; 5038965ea79SLois Curfman McInnes } 5048965ea79SLois Curfman McInnes 50539ddd567SLois Curfman McInnes static int MatView_MPIDense(PetscObject obj,Viewer viewer) 5068965ea79SLois Curfman McInnes { 5078965ea79SLois Curfman McInnes Mat mat = (Mat) obj; 50839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 5098965ea79SLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 51039ddd567SLois Curfman McInnes int ierr; 5118965ea79SLois Curfman McInnes 51239ddd567SLois Curfman McInnes if (!mdn->assembled) SETERRQ(1,"MatView_MPIDense:must assemble matrix"); 5138965ea79SLois Curfman McInnes if (!viewer) { 5148965ea79SLois Curfman McInnes viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 5158965ea79SLois Curfman McInnes } 51639ddd567SLois Curfman McInnes if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILE_VIEWER) { 51739ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 5188965ea79SLois Curfman McInnes } 5198965ea79SLois Curfman McInnes else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILES_VIEWER) { 52039ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 5218965ea79SLois Curfman McInnes } 5228965ea79SLois Curfman McInnes else if (vobj->type == BINARY_FILE_VIEWER) { 52339ddd567SLois Curfman McInnes return MatView_MPIDense_Binary(mat,viewer); 5248965ea79SLois Curfman McInnes } 5258965ea79SLois Curfman McInnes return 0; 5268965ea79SLois Curfman McInnes } 5278965ea79SLois Curfman McInnes 528*3501a2bdSLois Curfman McInnes static int MatGetInfo_MPIDense(Mat A,MatInfoType flag,int *nz, 5298965ea79SLois Curfman McInnes int *nzalloc,int *mem) 5308965ea79SLois Curfman McInnes { 531*3501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 532*3501a2bdSLois Curfman McInnes Mat mdn = mat->A; 53339ddd567SLois Curfman McInnes int ierr, isend[3], irecv[3]; 5348965ea79SLois Curfman McInnes 535*3501a2bdSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr); 5368965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 5378965ea79SLois Curfman McInnes *nz = isend[0]; *nzalloc = isend[1]; *mem = isend[2]; 5388965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 539*3501a2bdSLois Curfman McInnes MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_MAX,A->comm); 5408965ea79SLois Curfman McInnes *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2]; 5418965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 542*3501a2bdSLois Curfman McInnes MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_SUM,A->comm); 5438965ea79SLois Curfman McInnes *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2]; 5448965ea79SLois Curfman McInnes } 5458965ea79SLois Curfman McInnes return 0; 5468965ea79SLois Curfman McInnes } 5478965ea79SLois Curfman McInnes 54839ddd567SLois Curfman McInnes static int MatSetOption_MPIDense(Mat A,MatOption op) 5498965ea79SLois Curfman McInnes { 55039ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 5518965ea79SLois Curfman McInnes 5528965ea79SLois Curfman McInnes if (op == NO_NEW_NONZERO_LOCATIONS || 5538965ea79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 5548965ea79SLois Curfman McInnes op == COLUMNS_SORTED || 5558965ea79SLois Curfman McInnes op == ROW_ORIENTED) { 5568965ea79SLois Curfman McInnes MatSetOption(a->A,op); 5578965ea79SLois Curfman McInnes } 5588965ea79SLois Curfman McInnes else if (op == ROWS_SORTED || 5598965ea79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 5608965ea79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 5618965ea79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 56239ddd567SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_MPIDense:Option ignored\n"); 5638965ea79SLois Curfman McInnes else if (op == COLUMN_ORIENTED) 56439ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:COLUMN_ORIENTED");} 5658965ea79SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 56639ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:NO_NEW_DIAGONALS");} 5678965ea79SLois Curfman McInnes else 56839ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:unknown option");} 5698965ea79SLois Curfman McInnes return 0; 5708965ea79SLois Curfman McInnes } 5718965ea79SLois Curfman McInnes 572*3501a2bdSLois Curfman McInnes static int MatGetSize_MPIDense(Mat A,int *m,int *n) 5738965ea79SLois Curfman McInnes { 574*3501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5758965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 5768965ea79SLois Curfman McInnes return 0; 5778965ea79SLois Curfman McInnes } 5788965ea79SLois Curfman McInnes 579*3501a2bdSLois Curfman McInnes static int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 5808965ea79SLois Curfman McInnes { 581*3501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5828965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 5838965ea79SLois Curfman McInnes return 0; 5848965ea79SLois Curfman McInnes } 5858965ea79SLois Curfman McInnes 586*3501a2bdSLois Curfman McInnes static int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 5878965ea79SLois Curfman McInnes { 588*3501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5898965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 5908965ea79SLois Curfman McInnes return 0; 5918965ea79SLois Curfman McInnes } 5928965ea79SLois Curfman McInnes 593*3501a2bdSLois Curfman McInnes static int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 5948965ea79SLois Curfman McInnes { 595*3501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 59639ddd567SLois Curfman McInnes int lrow, rstart = mat->rstart, rend = mat->rend; 5978965ea79SLois Curfman McInnes 59839ddd567SLois Curfman McInnes if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIDense:only local rows") 5998965ea79SLois Curfman McInnes lrow = row - rstart; 60039ddd567SLois Curfman McInnes return MatGetRow(mat->A,lrow,nz,idx,v); 6018965ea79SLois Curfman McInnes } 6028965ea79SLois Curfman McInnes 60339ddd567SLois Curfman McInnes static int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 6048965ea79SLois Curfman McInnes { 6058965ea79SLois Curfman McInnes if (idx) PETSCFREE(*idx); 6068965ea79SLois Curfman McInnes if (v) PETSCFREE(*v); 6078965ea79SLois Curfman McInnes return 0; 6088965ea79SLois Curfman McInnes } 6098965ea79SLois Curfman McInnes 610*3501a2bdSLois Curfman McInnes static int MatNorm_MPIDense(Mat A,MatNormType type,double *norm) 611096963f5SLois Curfman McInnes { 612*3501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 613*3501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 614*3501a2bdSLois Curfman McInnes int ierr, i, j; 615*3501a2bdSLois Curfman McInnes double sum = 0.0; 616*3501a2bdSLois Curfman McInnes Scalar *v = mat->v; 617*3501a2bdSLois Curfman McInnes 618*3501a2bdSLois Curfman McInnes if (!mdn->assembled) SETERRQ(1,"MatNorm_MPIDense:Must assemble matrix"); 619*3501a2bdSLois Curfman McInnes if (mdn->size == 1) { 620*3501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm); CHKERRQ(ierr); 621*3501a2bdSLois Curfman McInnes } else { 622*3501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 623*3501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 624*3501a2bdSLois Curfman McInnes #if defined(PETSC_COMPLEX) 625*3501a2bdSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 626*3501a2bdSLois Curfman McInnes #else 627*3501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 628*3501a2bdSLois Curfman McInnes #endif 629*3501a2bdSLois Curfman McInnes } 630*3501a2bdSLois Curfman McInnes MPI_Allreduce((void*)&sum,(void*)norm,1,MPI_DOUBLE,MPI_SUM,A->comm); 631*3501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 632*3501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 633*3501a2bdSLois Curfman McInnes } 634*3501a2bdSLois Curfman McInnes else if (type == NORM_1) { 635*3501a2bdSLois Curfman McInnes double *tmp, *tmp2; 636*3501a2bdSLois Curfman McInnes tmp = (double *) PETSCMALLOC( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp); 637*3501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 638*3501a2bdSLois Curfman McInnes PetscZero(tmp,2*mdn->N*sizeof(double)); 639096963f5SLois Curfman McInnes *norm = 0.0; 640*3501a2bdSLois Curfman McInnes v = mat->v; 641*3501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 642*3501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 643*3501a2bdSLois Curfman McInnes #if defined(PETSC_COMPLEX) 644*3501a2bdSLois Curfman McInnes tmp[j] += abs(*v++); 645*3501a2bdSLois Curfman McInnes #else 646*3501a2bdSLois Curfman McInnes tmp[j] += fabs(*v++); 647*3501a2bdSLois Curfman McInnes #endif 648*3501a2bdSLois Curfman McInnes } 649*3501a2bdSLois Curfman McInnes } 650*3501a2bdSLois Curfman McInnes MPI_Allreduce((void*)tmp,(void*)tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm); 651*3501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 652*3501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 653*3501a2bdSLois Curfman McInnes } 654*3501a2bdSLois Curfman McInnes PETSCFREE(tmp); 655*3501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 656*3501a2bdSLois Curfman McInnes } 657*3501a2bdSLois Curfman McInnes else if (type == NORM_INFINITY) { /* max row norm */ 658*3501a2bdSLois Curfman McInnes double ntemp; 659*3501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr); 660*3501a2bdSLois Curfman McInnes MPI_Allreduce((void*)&ntemp,(void*)norm,1,MPI_DOUBLE,MPI_MAX,A->comm); 661*3501a2bdSLois Curfman McInnes } 662*3501a2bdSLois Curfman McInnes else { 663*3501a2bdSLois Curfman McInnes SETERRQ(1,"MatNorm_MPIDense:No support for two norm"); 664*3501a2bdSLois Curfman McInnes } 665*3501a2bdSLois Curfman McInnes } 666*3501a2bdSLois Curfman McInnes return 0; 667*3501a2bdSLois Curfman McInnes } 668*3501a2bdSLois Curfman McInnes 669*3501a2bdSLois Curfman McInnes static int MatTranspose_MPIDense(Mat A,Mat *matout) 670*3501a2bdSLois Curfman McInnes { 671*3501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 672*3501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 673*3501a2bdSLois Curfman McInnes Mat B; 674*3501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 675*3501a2bdSLois Curfman McInnes int j, i, ierr; 676*3501a2bdSLois Curfman McInnes Scalar *v; 677*3501a2bdSLois Curfman McInnes 678*3501a2bdSLois Curfman McInnes if (!matout && M != N) 679*3501a2bdSLois Curfman McInnes SETERRQ(1,"MatTranspose_MPIDense:Supports square matrix only in-place"); 680*3501a2bdSLois Curfman McInnes ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B); CHKERRQ(ierr); 681*3501a2bdSLois Curfman McInnes 682*3501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 683*3501a2bdSLois Curfman McInnes rwork = (int *) PETSCMALLOC(n*sizeof(int)); CHKPTRQ(rwork); 684*3501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 685*3501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 686*3501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr); 687*3501a2bdSLois Curfman McInnes v += m; 688*3501a2bdSLois Curfman McInnes } 689*3501a2bdSLois Curfman McInnes PETSCFREE(rwork); 690*3501a2bdSLois Curfman McInnes ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 691*3501a2bdSLois Curfman McInnes ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 692*3501a2bdSLois Curfman McInnes if (matout) { 693*3501a2bdSLois Curfman McInnes *matout = B; 694*3501a2bdSLois Curfman McInnes } else { 695*3501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 696*3501a2bdSLois Curfman McInnes PETSCFREE(a->rowners); 697*3501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 698*3501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 699*3501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 700*3501a2bdSLois Curfman McInnes PETSCFREE(a); 701*3501a2bdSLois Curfman McInnes PetscMemcpy(A,B,sizeof(struct _Mat)); 702*3501a2bdSLois Curfman McInnes PETSCHEADERDESTROY(B); 703*3501a2bdSLois Curfman McInnes } 704096963f5SLois Curfman McInnes return 0; 705096963f5SLois Curfman McInnes } 706096963f5SLois Curfman McInnes 70755659b69SBarry Smith static int MatCopyPrivate_MPIDense(Mat,Mat *,int); 7088965ea79SLois Curfman McInnes 7098965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 71039ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense, 71139ddd567SLois Curfman McInnes MatGetRow_MPIDense,MatRestoreRow_MPIDense, 71239ddd567SLois Curfman McInnes MatMult_MPIDense,MatMultAdd_MPIDense, 713096963f5SLois Curfman McInnes MatMultTrans_MPIDense,MatMultTransAdd_MPIDense, 71439ddd567SLois Curfman McInnes 0,0, 71539ddd567SLois Curfman McInnes 0,0,0, 716*3501a2bdSLois Curfman McInnes 0,0,MatTranspose_MPIDense, 71739ddd567SLois Curfman McInnes MatGetInfo_MPIDense,0, 718096963f5SLois Curfman McInnes MatGetDiagonal_MPIDense,0,MatNorm_MPIDense, 71939ddd567SLois Curfman McInnes MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense, 7208965ea79SLois Curfman McInnes 0, 72139ddd567SLois Curfman McInnes MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense, 72239ddd567SLois Curfman McInnes 0, 72339ddd567SLois Curfman McInnes 0,0,0,0, 72439ddd567SLois Curfman McInnes MatGetSize_MPIDense,MatGetLocalSize_MPIDense, 72539ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 72639ddd567SLois Curfman McInnes 0,0, 72739ddd567SLois Curfman McInnes 0,0,0,0,0,MatCopyPrivate_MPIDense}; 7288965ea79SLois Curfman McInnes 7298965ea79SLois Curfman McInnes /*@C 73039ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 7318965ea79SLois Curfman McInnes 7328965ea79SLois Curfman McInnes Input Parameters: 7338965ea79SLois Curfman McInnes . comm - MPI communicator 7348965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 7358965ea79SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated 7368965ea79SLois Curfman McInnes if N is given) 7378965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 7388965ea79SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated 7398965ea79SLois Curfman McInnes if n is given) 7408965ea79SLois Curfman McInnes 7418965ea79SLois Curfman McInnes Output Parameter: 7428965ea79SLois Curfman McInnes . newmat - the matrix 7438965ea79SLois Curfman McInnes 7448965ea79SLois Curfman McInnes Notes: 74539ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 74639ddd567SLois Curfman McInnes storage by columns. 7478965ea79SLois Curfman McInnes 7488965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 7498965ea79SLois Curfman McInnes (possibly both). 7508965ea79SLois Curfman McInnes 751*3501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 752*3501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 753*3501a2bdSLois Curfman McInnes 75439ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 7558965ea79SLois Curfman McInnes 75639ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 7578965ea79SLois Curfman McInnes @*/ 75839ddd567SLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Mat *newmat) 7598965ea79SLois Curfman McInnes { 7608965ea79SLois Curfman McInnes Mat mat; 76139ddd567SLois Curfman McInnes Mat_MPIDense *a; 76239ddd567SLois Curfman McInnes int ierr, i; 7638965ea79SLois Curfman McInnes 7648965ea79SLois Curfman McInnes *newmat = 0; 76539ddd567SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIDENSE,comm); 7668965ea79SLois Curfman McInnes PLogObjectCreate(mat); 76739ddd567SLois Curfman McInnes mat->data = (void *) (a = PETSCNEW(Mat_MPIDense)); CHKPTRQ(a); 7688965ea79SLois Curfman McInnes PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 76939ddd567SLois Curfman McInnes mat->destroy = MatDestroy_MPIDense; 77039ddd567SLois Curfman McInnes mat->view = MatView_MPIDense; 7718965ea79SLois Curfman McInnes mat->factor = 0; 7728965ea79SLois Curfman McInnes 7738965ea79SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 7748965ea79SLois Curfman McInnes MPI_Comm_rank(comm,&a->rank); 7758965ea79SLois Curfman McInnes MPI_Comm_size(comm,&a->size); 7768965ea79SLois Curfman McInnes 77739ddd567SLois Curfman McInnes if (M == PETSC_DECIDE) MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm); 7788965ea79SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);} 77939ddd567SLois Curfman McInnes 78039ddd567SLois Curfman McInnes /* each row stores all columns */ 78139ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 78239ddd567SLois Curfman McInnes if (n == PETSC_DECIDE) n = N; 78339ddd567SLois Curfman McInnes if (n != N) SETERRQ(1,"MatCreateMPIDense:For now, only n=N is supported"); 7848965ea79SLois Curfman McInnes a->N = N; 7858965ea79SLois Curfman McInnes a->M = M; 78639ddd567SLois Curfman McInnes a->m = m; 78739ddd567SLois Curfman McInnes a->n = n; 7888965ea79SLois Curfman McInnes 7898965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 79039ddd567SLois Curfman McInnes a->rowners = (int *) PETSCMALLOC((a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 79139ddd567SLois Curfman McInnes PLogObjectMemory(mat,(a->size+2)*sizeof(int)+sizeof(struct _Mat)+ 79239ddd567SLois Curfman McInnes sizeof(Mat_MPIDense)); 7938965ea79SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm); 7948965ea79SLois Curfman McInnes a->rowners[0] = 0; 7958965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 7968965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 7978965ea79SLois Curfman McInnes } 7988965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 7998965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 8008965ea79SLois Curfman McInnes 80139ddd567SLois Curfman McInnes ierr = MatCreateSeqDense(MPI_COMM_SELF,m,N,&a->A); CHKERRQ(ierr); 8028965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 8038965ea79SLois Curfman McInnes 8048965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 8058965ea79SLois Curfman McInnes ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr); 8068965ea79SLois Curfman McInnes 8078965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 8088965ea79SLois Curfman McInnes a->lvec = 0; 8098965ea79SLois Curfman McInnes a->Mvctx = 0; 8108965ea79SLois Curfman McInnes a->assembled = 0; 8118965ea79SLois Curfman McInnes 8128965ea79SLois Curfman McInnes *newmat = mat; 8138965ea79SLois Curfman McInnes return 0; 8148965ea79SLois Curfman McInnes } 8158965ea79SLois Curfman McInnes 816*3501a2bdSLois Curfman McInnes static int MatCopyPrivate_MPIDense(Mat A,Mat *newmat,int cpvalues) 8178965ea79SLois Curfman McInnes { 8188965ea79SLois Curfman McInnes Mat mat; 819*3501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 82039ddd567SLois Curfman McInnes int ierr; 8218965ea79SLois Curfman McInnes 82239ddd567SLois Curfman McInnes if (!oldmat->assembled) SETERRQ(1,"MatCopyPrivate_MPIDense:Must assemble matrix"); 8238965ea79SLois Curfman McInnes *newmat = 0; 824*3501a2bdSLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIDENSE,A->comm); 8258965ea79SLois Curfman McInnes PLogObjectCreate(mat); 82639ddd567SLois Curfman McInnes mat->data = (void *) (a = PETSCNEW(Mat_MPIDense)); CHKPTRQ(a); 8278965ea79SLois Curfman McInnes PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 82839ddd567SLois Curfman McInnes mat->destroy = MatDestroy_MPIDense; 82939ddd567SLois Curfman McInnes mat->view = MatView_MPIDense; 830*3501a2bdSLois Curfman McInnes mat->factor = A->factor; 8318965ea79SLois Curfman McInnes 8328965ea79SLois Curfman McInnes a->m = oldmat->m; 8338965ea79SLois Curfman McInnes a->n = oldmat->n; 8348965ea79SLois Curfman McInnes a->M = oldmat->M; 8358965ea79SLois Curfman McInnes a->N = oldmat->N; 8368965ea79SLois Curfman McInnes 8378965ea79SLois Curfman McInnes a->assembled = 1; 8388965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 8398965ea79SLois Curfman McInnes a->rend = oldmat->rend; 8408965ea79SLois Curfman McInnes a->size = oldmat->size; 8418965ea79SLois Curfman McInnes a->rank = oldmat->rank; 8428965ea79SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 8438965ea79SLois Curfman McInnes 8448965ea79SLois Curfman McInnes a->rowners = (int *) PETSCMALLOC((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 84539ddd567SLois Curfman McInnes PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense)); 8468965ea79SLois Curfman McInnes PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 8478965ea79SLois Curfman McInnes ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 8488965ea79SLois Curfman McInnes 8498965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 8508965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 85155659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 8528965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 8538965ea79SLois Curfman McInnes ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 8548965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 8558965ea79SLois Curfman McInnes *newmat = mat; 8568965ea79SLois Curfman McInnes return 0; 8578965ea79SLois Curfman McInnes } 8588965ea79SLois Curfman McInnes 8598965ea79SLois Curfman McInnes #include "sysio.h" 8608965ea79SLois Curfman McInnes 86139ddd567SLois Curfman McInnes int MatLoad_MPIDense(Viewer bview,MatType type,Mat *newmat) 8628965ea79SLois Curfman McInnes { 8638965ea79SLois Curfman McInnes Mat A; 8648965ea79SLois Curfman McInnes int i, nz, ierr, j,rstart, rend, fd; 8658965ea79SLois Curfman McInnes Scalar *vals,*svals; 8668965ea79SLois Curfman McInnes PetscObject vobj = (PetscObject) bview; 8678965ea79SLois Curfman McInnes MPI_Comm comm = vobj->comm; 8688965ea79SLois Curfman McInnes MPI_Status status; 8698965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 8708965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 8718965ea79SLois Curfman McInnes int tag = ((PetscObject)bview)->tag; 8728965ea79SLois Curfman McInnes 8738965ea79SLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 8748965ea79SLois Curfman McInnes if (!rank) { 8758965ea79SLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 8768965ea79SLois Curfman McInnes ierr = SYRead(fd,(char *)header,4,SYINT); CHKERRQ(ierr); 87739ddd567SLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:not matrix object"); 8788965ea79SLois Curfman McInnes } 8798965ea79SLois Curfman McInnes 8808965ea79SLois Curfman McInnes MPI_Bcast(header+1,3,MPI_INT,0,comm); 8818965ea79SLois Curfman McInnes M = header[1]; N = header[2]; 8828965ea79SLois Curfman McInnes /* determine ownership of all rows */ 8838965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 8848965ea79SLois Curfman McInnes rowners = (int *) PETSCMALLOC((size+2)*sizeof(int)); CHKPTRQ(rowners); 8858965ea79SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 8868965ea79SLois Curfman McInnes rowners[0] = 0; 8878965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 8888965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 8898965ea79SLois Curfman McInnes } 8908965ea79SLois Curfman McInnes rstart = rowners[rank]; 8918965ea79SLois Curfman McInnes rend = rowners[rank+1]; 8928965ea79SLois Curfman McInnes 8938965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 8948965ea79SLois Curfman McInnes ourlens = (int*) PETSCMALLOC( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 8958965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 8968965ea79SLois Curfman McInnes if (!rank) { 8978965ea79SLois Curfman McInnes rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths); 8988965ea79SLois Curfman McInnes ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 8998965ea79SLois Curfman McInnes sndcounts = (int*) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(sndcounts); 9008965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 9018965ea79SLois Curfman McInnes MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 9028965ea79SLois Curfman McInnes PETSCFREE(sndcounts); 9038965ea79SLois Curfman McInnes } 9048965ea79SLois Curfman McInnes else { 9058965ea79SLois Curfman McInnes MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 9068965ea79SLois Curfman McInnes } 9078965ea79SLois Curfman McInnes 9088965ea79SLois Curfman McInnes if (!rank) { 9098965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 9108965ea79SLois Curfman McInnes procsnz = (int*) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(procsnz); 9118965ea79SLois Curfman McInnes PetscZero(procsnz,size*sizeof(int)); 9128965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 9138965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 9148965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 9158965ea79SLois Curfman McInnes } 9168965ea79SLois Curfman McInnes } 9178965ea79SLois Curfman McInnes PETSCFREE(rowlengths); 9188965ea79SLois Curfman McInnes 9198965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 9208965ea79SLois Curfman McInnes maxnz = 0; 9218965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 9228965ea79SLois Curfman McInnes maxnz = PETSCMAX(maxnz,procsnz[i]); 9238965ea79SLois Curfman McInnes } 9248965ea79SLois Curfman McInnes cols = (int *) PETSCMALLOC( maxnz*sizeof(int) ); CHKPTRQ(cols); 9258965ea79SLois Curfman McInnes 9268965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 9278965ea79SLois Curfman McInnes nz = procsnz[0]; 9288965ea79SLois Curfman McInnes mycols = (int *) PETSCMALLOC( nz*sizeof(int) ); CHKPTRQ(mycols); 9298965ea79SLois Curfman McInnes ierr = SYRead(fd,mycols,nz,SYINT); CHKERRQ(ierr); 9308965ea79SLois Curfman McInnes 9318965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 9328965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 9338965ea79SLois Curfman McInnes nz = procsnz[i]; 9348965ea79SLois Curfman McInnes ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr); 9358965ea79SLois Curfman McInnes MPI_Send(cols,nz,MPI_INT,i,tag,comm); 9368965ea79SLois Curfman McInnes } 9378965ea79SLois Curfman McInnes PETSCFREE(cols); 9388965ea79SLois Curfman McInnes } 9398965ea79SLois Curfman McInnes else { 9408965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 9418965ea79SLois Curfman McInnes nz = 0; 9428965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 9438965ea79SLois Curfman McInnes nz += ourlens[i]; 9448965ea79SLois Curfman McInnes } 9458965ea79SLois Curfman McInnes mycols = (int*) PETSCMALLOC( nz*sizeof(int) ); CHKPTRQ(mycols); 9468965ea79SLois Curfman McInnes 9478965ea79SLois Curfman McInnes /* receive message of column indices*/ 9488965ea79SLois Curfman McInnes MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 9498965ea79SLois Curfman McInnes MPI_Get_count(&status,MPI_INT,&maxnz); 95039ddd567SLois Curfman McInnes if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file"); 9518965ea79SLois Curfman McInnes } 9528965ea79SLois Curfman McInnes 9538965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 9548965ea79SLois Curfman McInnes PetscZero(offlens,m*sizeof(int)); 9558965ea79SLois Curfman McInnes jj = 0; 9568965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 9578965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 9588965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 9598965ea79SLois Curfman McInnes jj++; 9608965ea79SLois Curfman McInnes } 9618965ea79SLois Curfman McInnes } 9628965ea79SLois Curfman McInnes 9638965ea79SLois Curfman McInnes /* create our matrix */ 9648965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 9658965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 9668965ea79SLois Curfman McInnes } 96739ddd567SLois Curfman McInnes if (type == MATMPIDENSE) { 96839ddd567SLois Curfman McInnes ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 9698965ea79SLois Curfman McInnes } 9708965ea79SLois Curfman McInnes A = *newmat; 9718965ea79SLois Curfman McInnes MatSetOption(A,COLUMNS_SORTED); 9728965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 9738965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 9748965ea79SLois Curfman McInnes } 9758965ea79SLois Curfman McInnes 9768965ea79SLois Curfman McInnes if (!rank) { 9778965ea79SLois Curfman McInnes vals = (Scalar *) PETSCMALLOC( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 9788965ea79SLois Curfman McInnes 9798965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 9808965ea79SLois Curfman McInnes nz = procsnz[0]; 9818965ea79SLois Curfman McInnes ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 9828965ea79SLois Curfman McInnes 9838965ea79SLois Curfman McInnes /* insert into matrix */ 9848965ea79SLois Curfman McInnes jj = rstart; 9858965ea79SLois Curfman McInnes smycols = mycols; 9868965ea79SLois Curfman McInnes svals = vals; 9878965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 9888965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 9898965ea79SLois Curfman McInnes smycols += ourlens[i]; 9908965ea79SLois Curfman McInnes svals += ourlens[i]; 9918965ea79SLois Curfman McInnes jj++; 9928965ea79SLois Curfman McInnes } 9938965ea79SLois Curfman McInnes 9948965ea79SLois Curfman McInnes /* read in other processors and ship out */ 9958965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 9968965ea79SLois Curfman McInnes nz = procsnz[i]; 9978965ea79SLois Curfman McInnes ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 9988965ea79SLois Curfman McInnes MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 9998965ea79SLois Curfman McInnes } 10008965ea79SLois Curfman McInnes PETSCFREE(procsnz); 10018965ea79SLois Curfman McInnes } 10028965ea79SLois Curfman McInnes else { 10038965ea79SLois Curfman McInnes /* receive numeric values */ 10048965ea79SLois Curfman McInnes vals = (Scalar*) PETSCMALLOC( nz*sizeof(Scalar) ); CHKPTRQ(vals); 10058965ea79SLois Curfman McInnes 10068965ea79SLois Curfman McInnes /* receive message of values*/ 10078965ea79SLois Curfman McInnes MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 10088965ea79SLois Curfman McInnes MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 100939ddd567SLois Curfman McInnes if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file"); 10108965ea79SLois Curfman McInnes 10118965ea79SLois Curfman McInnes /* insert into matrix */ 10128965ea79SLois Curfman McInnes jj = rstart; 10138965ea79SLois Curfman McInnes smycols = mycols; 10148965ea79SLois Curfman McInnes svals = vals; 10158965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 10168965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 10178965ea79SLois Curfman McInnes smycols += ourlens[i]; 10188965ea79SLois Curfman McInnes svals += ourlens[i]; 10198965ea79SLois Curfman McInnes jj++; 10208965ea79SLois Curfman McInnes } 10218965ea79SLois Curfman McInnes } 10228965ea79SLois Curfman McInnes PETSCFREE(ourlens); PETSCFREE(vals); PETSCFREE(mycols); PETSCFREE(rowners); 10238965ea79SLois Curfman McInnes 10248965ea79SLois Curfman McInnes ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 10258965ea79SLois Curfman McInnes ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 10268965ea79SLois Curfman McInnes return 0; 10278965ea79SLois Curfman McInnes } 1028