1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*549d3d68SSatish Balay static char vcid[] = "$Id: mpidense.c,v 1.112 1999/04/19 22:11:55 bsmith Exp balay $"; 38965ea79SLois Curfman McInnes #endif 48965ea79SLois Curfman McInnes 5ed3cc1f0SBarry Smith /* 6ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 7ed3cc1f0SBarry Smith */ 8ed3cc1f0SBarry Smith 970f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 10f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 118965ea79SLois Curfman McInnes 127ef1d9bdSSatish Balay extern int MatSetUpMultiply_MPIDense(Mat); 137ef1d9bdSSatish Balay 145615d1e5SSatish Balay #undef __FUNC__ 155615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIDense" 168f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 178965ea79SLois Curfman McInnes { 1839b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense *) mat->data; 1939b7565bSBarry Smith int ierr, i, j, rstart = A->rstart, rend = A->rend, row; 2039b7565bSBarry Smith int roworiented = A->roworiented; 218965ea79SLois Curfman McInnes 223a40ed3dSBarry Smith PetscFunctionBegin; 238965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 245ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 25a8c6a408SBarry Smith if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 268965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 278965ea79SLois Curfman McInnes row = idxm[i] - rstart; 2839b7565bSBarry Smith if (roworiented) { 2939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 303a40ed3dSBarry Smith } else { 318965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 325ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 33a8c6a408SBarry Smith if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 3439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 3539b7565bSBarry Smith } 368965ea79SLois Curfman McInnes } 373a40ed3dSBarry Smith } else { 383782ba37SSatish Balay if ( !A->donotstash) { 3939b7565bSBarry Smith if (roworiented) { 408798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 41d36fbae8SSatish Balay } else { 428798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 4339b7565bSBarry Smith } 44b49de8d1SLois Curfman McInnes } 45b49de8d1SLois Curfman McInnes } 463782ba37SSatish Balay } 473a40ed3dSBarry Smith PetscFunctionReturn(0); 48b49de8d1SLois Curfman McInnes } 49b49de8d1SLois Curfman McInnes 505615d1e5SSatish Balay #undef __FUNC__ 515615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIDense" 528f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 53b49de8d1SLois Curfman McInnes { 54b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 55b49de8d1SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 56b49de8d1SLois Curfman McInnes 573a40ed3dSBarry Smith PetscFunctionBegin; 58b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 59a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 60a8c6a408SBarry Smith if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 61b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 62b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 63b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 64a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 65a8c6a408SBarry Smith if (idxn[j] >= mdn->N) { 66a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 67a8c6a408SBarry Smith } 68b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 69b49de8d1SLois Curfman McInnes } 70a8c6a408SBarry Smith } else { 71a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 728965ea79SLois Curfman McInnes } 738965ea79SLois Curfman McInnes } 743a40ed3dSBarry Smith PetscFunctionReturn(0); 758965ea79SLois Curfman McInnes } 768965ea79SLois Curfman McInnes 775615d1e5SSatish Balay #undef __FUNC__ 785615d1e5SSatish Balay #define __FUNC__ "MatGetArray_MPIDense" 798f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array) 80ff14e315SSatish Balay { 81ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense *) A->data; 82ff14e315SSatish Balay int ierr; 83ff14e315SSatish Balay 843a40ed3dSBarry Smith PetscFunctionBegin; 85ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 863a40ed3dSBarry Smith PetscFunctionReturn(0); 87ff14e315SSatish Balay } 88ff14e315SSatish Balay 895615d1e5SSatish Balay #undef __FUNC__ 905615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_MPIDense" 918f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array) 92ff14e315SSatish Balay { 933a40ed3dSBarry Smith PetscFunctionBegin; 943a40ed3dSBarry Smith PetscFunctionReturn(0); 95ff14e315SSatish Balay } 96ff14e315SSatish Balay 975615d1e5SSatish Balay #undef __FUNC__ 985615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIDense" 998f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1008965ea79SLois Curfman McInnes { 10139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1028965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 103d36fbae8SSatish Balay int ierr,nstash,reallocs; 1048965ea79SLois Curfman McInnes InsertMode addv; 1058965ea79SLois Curfman McInnes 1063a40ed3dSBarry Smith PetscFunctionBegin; 1078965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 108ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1097056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 110a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs"); 1118965ea79SLois Curfman McInnes } 112e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1138965ea79SLois Curfman McInnes 1148798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1158798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 116d36fbae8SSatish Balay PLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n", 117d36fbae8SSatish Balay nstash,reallocs); 1183a40ed3dSBarry Smith PetscFunctionReturn(0); 1198965ea79SLois Curfman McInnes } 1208965ea79SLois Curfman McInnes 1215615d1e5SSatish Balay #undef __FUNC__ 1225615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense" 1238f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1248965ea79SLois Curfman McInnes { 12539ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 1267ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 1277ef1d9bdSSatish Balay Scalar *val; 128e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 1298965ea79SLois Curfman McInnes 1303a40ed3dSBarry Smith PetscFunctionBegin; 1318965ea79SLois Curfman McInnes /* wait on receives */ 1327ef1d9bdSSatish Balay while (1) { 1338798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 1347ef1d9bdSSatish Balay if (!flg) break; 1358965ea79SLois Curfman McInnes 1367ef1d9bdSSatish Balay for ( i=0; i<n; ) { 1377ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 1387ef1d9bdSSatish Balay for ( j=i,rstart=row[j]; j<n; j++ ) { if (row[j] != rstart) break; } 1397ef1d9bdSSatish Balay if (j < n) ncols = j-i; 1407ef1d9bdSSatish Balay else ncols = n-i; 1417ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 1427ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 1437ef1d9bdSSatish Balay i = j; 1448965ea79SLois Curfman McInnes } 1457ef1d9bdSSatish Balay } 1468798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 1478965ea79SLois Curfman McInnes 14839ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 14939ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 1508965ea79SLois Curfman McInnes 1516d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 15239ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 1538965ea79SLois Curfman McInnes } 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 1558965ea79SLois Curfman McInnes } 1568965ea79SLois Curfman McInnes 1575615d1e5SSatish Balay #undef __FUNC__ 1585615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense" 1598f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 1608965ea79SLois Curfman McInnes { 1613a40ed3dSBarry Smith int ierr; 16239ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 1633a40ed3dSBarry Smith 1643a40ed3dSBarry Smith PetscFunctionBegin; 1653a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 1663a40ed3dSBarry Smith PetscFunctionReturn(0); 1678965ea79SLois Curfman McInnes } 1688965ea79SLois Curfman McInnes 1695615d1e5SSatish Balay #undef __FUNC__ 1705615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense" 1718f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 1724e220ebcSLois Curfman McInnes { 1733a40ed3dSBarry Smith PetscFunctionBegin; 1744e220ebcSLois Curfman McInnes *bs = 1; 1753a40ed3dSBarry Smith PetscFunctionReturn(0); 1764e220ebcSLois Curfman McInnes } 1774e220ebcSLois Curfman McInnes 1788965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 1798965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 1808965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 1813501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 1828965ea79SLois Curfman McInnes routine. 1838965ea79SLois Curfman McInnes */ 1845615d1e5SSatish Balay #undef __FUNC__ 1855615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense" 1868f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 1878965ea79SLois Curfman McInnes { 18839ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 1898965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 1908965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 1918965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 1928965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 1938965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 1948965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 1958965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 1968965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 1978965ea79SLois Curfman McInnes IS istmp; 1988965ea79SLois Curfman McInnes 1993a40ed3dSBarry Smith PetscFunctionBegin; 20077c4ece6SBarry Smith ierr = ISGetSize(is,&N);CHKERRQ(ierr); 2018965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2028965ea79SLois Curfman McInnes 2038965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2040452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(nprocs); 205*549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 206*549d3d68SSatish Balay procs = nprocs + size; 2070452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/ 2088965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2098965ea79SLois Curfman McInnes idx = rows[i]; 2108965ea79SLois Curfman McInnes found = 0; 2118965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2128965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2138965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2148965ea79SLois Curfman McInnes } 2158965ea79SLois Curfman McInnes } 216a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 2178965ea79SLois Curfman McInnes } 2188965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 2198965ea79SLois Curfman McInnes 2208965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 2210452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) );CHKPTRQ(work); 222ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2238965ea79SLois Curfman McInnes nrecvs = work[rank]; 224ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 2258965ea79SLois Curfman McInnes nmax = work[rank]; 2260452661fSBarry Smith PetscFree(work); 2278965ea79SLois Curfman McInnes 2288965ea79SLois Curfman McInnes /* post receives: */ 2293a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 2303a40ed3dSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 2318965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 232ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 2338965ea79SLois Curfman McInnes } 2348965ea79SLois Curfman McInnes 2358965ea79SLois Curfman McInnes /* do sends: 2368965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 2378965ea79SLois Curfman McInnes the ith processor 2388965ea79SLois Curfman McInnes */ 2390452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(svalues); 2407056b6fcSBarry Smith send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 2410452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) );CHKPTRQ(starts); 2428965ea79SLois Curfman McInnes starts[0] = 0; 2438965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2448965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2458965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 2468965ea79SLois Curfman McInnes } 2478965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 2488965ea79SLois Curfman McInnes 2498965ea79SLois Curfman McInnes starts[0] = 0; 2508965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2518965ea79SLois Curfman McInnes count = 0; 2528965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 2538965ea79SLois Curfman McInnes if (procs[i]) { 254ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 2558965ea79SLois Curfman McInnes } 2568965ea79SLois Curfman McInnes } 2570452661fSBarry Smith PetscFree(starts); 2588965ea79SLois Curfman McInnes 2598965ea79SLois Curfman McInnes base = owners[rank]; 2608965ea79SLois Curfman McInnes 2618965ea79SLois Curfman McInnes /* wait on receives */ 2620452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) );CHKPTRQ(lens); 2638965ea79SLois Curfman McInnes source = lens + nrecvs; 2648965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 2658965ea79SLois Curfman McInnes while (count) { 266ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 2678965ea79SLois Curfman McInnes /* unpack receives into our local space */ 268ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 2698965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 2708965ea79SLois Curfman McInnes lens[imdex] = n; 2718965ea79SLois Curfman McInnes slen += n; 2728965ea79SLois Curfman McInnes count--; 2738965ea79SLois Curfman McInnes } 2740452661fSBarry Smith PetscFree(recv_waits); 2758965ea79SLois Curfman McInnes 2768965ea79SLois Curfman McInnes /* move the data into the send scatter */ 2770452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) );CHKPTRQ(lrows); 2788965ea79SLois Curfman McInnes count = 0; 2798965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2808965ea79SLois Curfman McInnes values = rvalues + i*nmax; 2818965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 2828965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 2838965ea79SLois Curfman McInnes } 2848965ea79SLois Curfman McInnes } 2850452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 2860452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 2878965ea79SLois Curfman McInnes 2888965ea79SLois Curfman McInnes /* actually zap the local rows */ 289029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 2908965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 2910452661fSBarry Smith PetscFree(lrows); 2928965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 2938965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 2948965ea79SLois Curfman McInnes 2958965ea79SLois Curfman McInnes /* wait on sends */ 2968965ea79SLois Curfman McInnes if (nsends) { 2977056b6fcSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 298ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 2990452661fSBarry Smith PetscFree(send_status); 3008965ea79SLois Curfman McInnes } 3010452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 3028965ea79SLois Curfman McInnes 3033a40ed3dSBarry Smith PetscFunctionReturn(0); 3048965ea79SLois Curfman McInnes } 3058965ea79SLois Curfman McInnes 3065615d1e5SSatish Balay #undef __FUNC__ 3075615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense" 3088f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3098965ea79SLois Curfman McInnes { 31039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3118965ea79SLois Curfman McInnes int ierr; 312c456f294SBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 31443a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 31543a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 31644cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3173a40ed3dSBarry Smith PetscFunctionReturn(0); 3188965ea79SLois Curfman McInnes } 3198965ea79SLois Curfman McInnes 3205615d1e5SSatish Balay #undef __FUNC__ 3215615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense" 3228f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3238965ea79SLois Curfman McInnes { 32439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3258965ea79SLois Curfman McInnes int ierr; 326c456f294SBarry Smith 3273a40ed3dSBarry Smith PetscFunctionBegin; 32843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 32943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 33044cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 3313a40ed3dSBarry Smith PetscFunctionReturn(0); 3328965ea79SLois Curfman McInnes } 3338965ea79SLois Curfman McInnes 3345615d1e5SSatish Balay #undef __FUNC__ 3355615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense" 3368f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 337096963f5SLois Curfman McInnes { 338096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 339096963f5SLois Curfman McInnes int ierr; 3403501a2bdSLois Curfman McInnes Scalar zero = 0.0; 341096963f5SLois Curfman McInnes 3423a40ed3dSBarry Smith PetscFunctionBegin; 3433501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 34444cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 345537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 346537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 3473a40ed3dSBarry Smith PetscFunctionReturn(0); 348096963f5SLois Curfman McInnes } 349096963f5SLois Curfman McInnes 3505615d1e5SSatish Balay #undef __FUNC__ 3515615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense" 3528f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 353096963f5SLois Curfman McInnes { 354096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 355096963f5SLois Curfman McInnes int ierr; 356096963f5SLois Curfman McInnes 3573a40ed3dSBarry Smith PetscFunctionBegin; 3583501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 35944cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 360537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 361537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 3623a40ed3dSBarry Smith PetscFunctionReturn(0); 363096963f5SLois Curfman McInnes } 364096963f5SLois Curfman McInnes 3655615d1e5SSatish Balay #undef __FUNC__ 3665615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense" 3678f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 3688965ea79SLois Curfman McInnes { 36939ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 370096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 37144cd7ae7SLois Curfman McInnes int ierr, len, i, n, m = a->m, radd; 37244cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 373ed3cc1f0SBarry Smith 3743a40ed3dSBarry Smith PetscFunctionBegin; 37544cd7ae7SLois Curfman McInnes VecSet(&zero,v); 376096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x);CHKERRQ(ierr); 377096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 378a8c6a408SBarry Smith if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 37944cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 3807ddc982cSLois Curfman McInnes radd = a->rstart*m; 38144cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 382096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 383096963f5SLois Curfman McInnes } 3843a40ed3dSBarry Smith PetscFunctionReturn(0); 3858965ea79SLois Curfman McInnes } 3868965ea79SLois Curfman McInnes 3875615d1e5SSatish Balay #undef __FUNC__ 3885615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense" 389e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 3908965ea79SLois Curfman McInnes { 3913501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3928965ea79SLois Curfman McInnes int ierr; 393ed3cc1f0SBarry Smith 3943a40ed3dSBarry Smith PetscFunctionBegin; 39594d884c6SBarry Smith if (--mat->refct > 0) PetscFunctionReturn(0); 39694d884c6SBarry Smith 39794d884c6SBarry Smith if (mat->mapping) { 39894d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 39994d884c6SBarry Smith } 40094d884c6SBarry Smith if (mat->bmapping) { 40194d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 40294d884c6SBarry Smith } 4033a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 404e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4058965ea79SLois Curfman McInnes #endif 4068798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 4070452661fSBarry Smith PetscFree(mdn->rowners); 4083501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4093501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4103501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 411622d7880SLois Curfman McInnes if (mdn->factor) { 412622d7880SLois Curfman McInnes if (mdn->factor->temp) PetscFree(mdn->factor->temp); 413622d7880SLois Curfman McInnes if (mdn->factor->tag) PetscFree(mdn->factor->tag); 414622d7880SLois Curfman McInnes if (mdn->factor->pivots) PetscFree(mdn->factor->pivots); 415622d7880SLois Curfman McInnes PetscFree(mdn->factor); 416622d7880SLois Curfman McInnes } 4170452661fSBarry Smith PetscFree(mdn); 41861b13de0SBarry Smith if (mat->rmap) { 41961b13de0SBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 42061b13de0SBarry Smith } 42161b13de0SBarry Smith if (mat->cmap) { 42261b13de0SBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 42390f02eecSBarry Smith } 4248965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4250452661fSBarry Smith PetscHeaderDestroy(mat); 4263a40ed3dSBarry Smith PetscFunctionReturn(0); 4278965ea79SLois Curfman McInnes } 42839ddd567SLois Curfman McInnes 4295615d1e5SSatish Balay #undef __FUNC__ 4305615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary" 43139ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4328965ea79SLois Curfman McInnes { 43339ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4348965ea79SLois Curfman McInnes int ierr; 4357056b6fcSBarry Smith 4363a40ed3dSBarry Smith PetscFunctionBegin; 43739ddd567SLois Curfman McInnes if (mdn->size == 1) { 43839ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4398965ea79SLois Curfman McInnes } 440a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 4413a40ed3dSBarry Smith PetscFunctionReturn(0); 4428965ea79SLois Curfman McInnes } 4438965ea79SLois Curfman McInnes 4445615d1e5SSatish Balay #undef __FUNC__ 4455615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII" 44639ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 4478965ea79SLois Curfman McInnes { 44839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 44977ed5343SBarry Smith int ierr, format, size = mdn->size, rank = mdn->rank; 4508965ea79SLois Curfman McInnes FILE *fd; 45119bcc07fSBarry Smith ViewerType vtype; 4528965ea79SLois Curfman McInnes 4533a40ed3dSBarry Smith PetscFunctionBegin; 4543a40ed3dSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 45590ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 45690ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 457639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 4584e220ebcSLois Curfman McInnes MatInfo info; 4594e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 46077c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 4614e220ebcSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 4624e220ebcSLois Curfman McInnes (int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 463096963f5SLois Curfman McInnes fflush(fd); 46477c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 4653501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 4663a40ed3dSBarry Smith PetscFunctionReturn(0); 46796f6c058SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 4698965ea79SLois Curfman McInnes } 47077ed5343SBarry Smith 4718965ea79SLois Curfman McInnes if (size == 1) { 47239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4733a40ed3dSBarry Smith } else { 4748965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 4758965ea79SLois Curfman McInnes Mat A; 47639ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 47739ddd567SLois Curfman McInnes Scalar *vals; 47839ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 4798965ea79SLois Curfman McInnes 4808965ea79SLois Curfman McInnes if (!rank) { 4810513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 4823a40ed3dSBarry Smith } else { 4830513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 4848965ea79SLois Curfman McInnes } 4858965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 4868965ea79SLois Curfman McInnes 48739ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 48839ddd567SLois Curfman McInnes but it's quick for now */ 48939ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 4908965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 49139ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 49239ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 49339ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 49439ddd567SLois Curfman McInnes row++; 4958965ea79SLois Curfman McInnes } 4968965ea79SLois Curfman McInnes 4976d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4986d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4998965ea79SLois Curfman McInnes if (!rank) { 50039ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer);CHKERRQ(ierr); 5018965ea79SLois Curfman McInnes } 5028965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5038965ea79SLois Curfman McInnes } 5043a40ed3dSBarry Smith PetscFunctionReturn(0); 5058965ea79SLois Curfman McInnes } 5068965ea79SLois Curfman McInnes 5075615d1e5SSatish Balay #undef __FUNC__ 5085615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense" 509e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer) 5108965ea79SLois Curfman McInnes { 51139ddd567SLois Curfman McInnes int ierr; 512bcd2baecSBarry Smith ViewerType vtype; 5138965ea79SLois Curfman McInnes 514bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5153f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 51639ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer);CHKERRQ(ierr); 5173f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5183a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 5195cd90555SBarry Smith } else { 5205cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5218965ea79SLois Curfman McInnes } 5223a40ed3dSBarry Smith PetscFunctionReturn(0); 5238965ea79SLois Curfman McInnes } 5248965ea79SLois Curfman McInnes 5255615d1e5SSatish Balay #undef __FUNC__ 5265615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense" 5278f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 5288965ea79SLois Curfman McInnes { 5293501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5303501a2bdSLois Curfman McInnes Mat mdn = mat->A; 5314e220ebcSLois Curfman McInnes int ierr; 5324e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 5338965ea79SLois Curfman McInnes 5343a40ed3dSBarry Smith PetscFunctionBegin; 5354e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 5364e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 5374e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 5384e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 5394e220ebcSLois Curfman McInnes info->block_size = 1.0; 5404e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 5414e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 5424e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 5438965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 5444e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 5454e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 5464e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 5474e220ebcSLois Curfman McInnes info->memory = isend[3]; 5484e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 5498965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 550f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 5514e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 5524e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 5534e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 5544e220ebcSLois Curfman McInnes info->memory = irecv[3]; 5554e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 5568965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 557f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 5584e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 5594e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 5604e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 5614e220ebcSLois Curfman McInnes info->memory = irecv[3]; 5624e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 5638965ea79SLois Curfman McInnes } 5644e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 5654e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 5664e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 5673a40ed3dSBarry Smith PetscFunctionReturn(0); 5688965ea79SLois Curfman McInnes } 5698965ea79SLois Curfman McInnes 5708c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 5718aaee692SLois Curfman McInnes extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 5728aaee692SLois Curfman McInnes extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 5738aaee692SLois Curfman McInnes extern int MatSolve_MPIDense(Mat,Vec,Vec); 5748c469469SLois Curfman McInnes extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 5758aaee692SLois Curfman McInnes extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 5768aaee692SLois Curfman McInnes extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 5778aaee692SLois Curfman McInnes 5785615d1e5SSatish Balay #undef __FUNC__ 5795615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense" 5808f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 5818965ea79SLois Curfman McInnes { 58239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 5838965ea79SLois Curfman McInnes 5843a40ed3dSBarry Smith PetscFunctionBegin; 5856d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 5866d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 5874787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 5884787f768SSatish Balay op == MAT_NEW_NONZERO_ALLOCATION_ERR || 589219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 590219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED) { 591b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 592b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 593aeafbbfcSLois Curfman McInnes a->roworiented = 1; 5948965ea79SLois Curfman McInnes MatSetOption(a->A,op); 595b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 596219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 5976d4a8577SBarry Smith op == MAT_SYMMETRIC || 5986d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 599b51ba29fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 600b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 601981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 6023a40ed3dSBarry Smith } else if (op == MAT_COLUMN_ORIENTED) { 6033a40ed3dSBarry Smith a->roworiented = 0; MatSetOption(a->A,op); 6043782ba37SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 6053782ba37SSatish Balay a->donotstash = 1; 6063a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) { 6073a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 6083a40ed3dSBarry Smith } else { 6093a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 6103a40ed3dSBarry Smith } 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6128965ea79SLois Curfman McInnes } 6138965ea79SLois Curfman McInnes 6145615d1e5SSatish Balay #undef __FUNC__ 6155615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense" 6168f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n) 6178965ea79SLois Curfman McInnes { 6183501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6193a40ed3dSBarry Smith 6203a40ed3dSBarry Smith PetscFunctionBegin; 6218965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 6223a40ed3dSBarry Smith PetscFunctionReturn(0); 6238965ea79SLois Curfman McInnes } 6248965ea79SLois Curfman McInnes 6255615d1e5SSatish Balay #undef __FUNC__ 6265615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense" 6278f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 6288965ea79SLois Curfman McInnes { 6293501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6303a40ed3dSBarry Smith 6313a40ed3dSBarry Smith PetscFunctionBegin; 6328965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 6333a40ed3dSBarry Smith PetscFunctionReturn(0); 6348965ea79SLois Curfman McInnes } 6358965ea79SLois Curfman McInnes 6365615d1e5SSatish Balay #undef __FUNC__ 6375615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense" 6388f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 6398965ea79SLois Curfman McInnes { 6403501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6413a40ed3dSBarry Smith 6423a40ed3dSBarry Smith PetscFunctionBegin; 6438965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 6443a40ed3dSBarry Smith PetscFunctionReturn(0); 6458965ea79SLois Curfman McInnes } 6468965ea79SLois Curfman McInnes 6475615d1e5SSatish Balay #undef __FUNC__ 6485615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense" 6498f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 6508965ea79SLois Curfman McInnes { 6513501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6523a40ed3dSBarry Smith int lrow, rstart = mat->rstart, rend = mat->rend,ierr; 6538965ea79SLois Curfman McInnes 6543a40ed3dSBarry Smith PetscFunctionBegin; 655a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 6568965ea79SLois Curfman McInnes lrow = row - rstart; 6573a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 6583a40ed3dSBarry Smith PetscFunctionReturn(0); 6598965ea79SLois Curfman McInnes } 6608965ea79SLois Curfman McInnes 6615615d1e5SSatish Balay #undef __FUNC__ 6625615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense" 6638f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 6648965ea79SLois Curfman McInnes { 6653a40ed3dSBarry Smith PetscFunctionBegin; 6660452661fSBarry Smith if (idx) PetscFree(*idx); 6670452661fSBarry Smith if (v) PetscFree(*v); 6683a40ed3dSBarry Smith PetscFunctionReturn(0); 6698965ea79SLois Curfman McInnes } 6708965ea79SLois Curfman McInnes 6715615d1e5SSatish Balay #undef __FUNC__ 6725615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense" 6738f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm) 674096963f5SLois Curfman McInnes { 6753501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 6763501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 6773501a2bdSLois Curfman McInnes int ierr, i, j; 6783501a2bdSLois Curfman McInnes double sum = 0.0; 6793501a2bdSLois Curfman McInnes Scalar *v = mat->v; 6803501a2bdSLois Curfman McInnes 6813a40ed3dSBarry Smith PetscFunctionBegin; 6823501a2bdSLois Curfman McInnes if (mdn->size == 1) { 6833501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm);CHKERRQ(ierr); 6843501a2bdSLois Curfman McInnes } else { 6853501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 6863501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 6873a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 688e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 6893501a2bdSLois Curfman McInnes #else 6903501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 6913501a2bdSLois Curfman McInnes #endif 6923501a2bdSLois Curfman McInnes } 693ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 6943501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 6953501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 6963a40ed3dSBarry Smith } else if (type == NORM_1) { 6973501a2bdSLois Curfman McInnes double *tmp, *tmp2; 6980452661fSBarry Smith tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) );CHKPTRQ(tmp); 6993501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 700*549d3d68SSatish Balay ierr = PetscMemzero(tmp,2*mdn->N*sizeof(double));CHKERRQ(ierr); 701096963f5SLois Curfman McInnes *norm = 0.0; 7023501a2bdSLois Curfman McInnes v = mat->v; 7033501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 7043501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 70567e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7063501a2bdSLois Curfman McInnes } 7073501a2bdSLois Curfman McInnes } 708ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7093501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 7103501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 7113501a2bdSLois Curfman McInnes } 7120452661fSBarry Smith PetscFree(tmp); 7133501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 7143a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 7153501a2bdSLois Curfman McInnes double ntemp; 7163501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 717ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 7183a40ed3dSBarry Smith } else { 719a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 7203501a2bdSLois Curfman McInnes } 7213501a2bdSLois Curfman McInnes } 7223a40ed3dSBarry Smith PetscFunctionReturn(0); 7233501a2bdSLois Curfman McInnes } 7243501a2bdSLois Curfman McInnes 7255615d1e5SSatish Balay #undef __FUNC__ 7265615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense" 7278f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 7283501a2bdSLois Curfman McInnes { 7293501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 7303501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 7313501a2bdSLois Curfman McInnes Mat B; 7323501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 7333501a2bdSLois Curfman McInnes int j, i, ierr; 7343501a2bdSLois Curfman McInnes Scalar *v; 7353501a2bdSLois Curfman McInnes 7363a40ed3dSBarry Smith PetscFunctionBegin; 7377056b6fcSBarry Smith if (matout == PETSC_NULL && M != N) { 738a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 7397056b6fcSBarry Smith } 7407056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 7413501a2bdSLois Curfman McInnes 7423501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 7430452661fSBarry Smith rwork = (int *) PetscMalloc(n*sizeof(int));CHKPTRQ(rwork); 7443501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 7453501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 7463501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 7473501a2bdSLois Curfman McInnes v += m; 7483501a2bdSLois Curfman McInnes } 7490452661fSBarry Smith PetscFree(rwork); 7506d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7516d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7523638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 7533501a2bdSLois Curfman McInnes *matout = B; 7543501a2bdSLois Curfman McInnes } else { 755f830108cSBarry Smith PetscOps *Abops; 75609dc0095SBarry Smith MatOps Aops; 757f830108cSBarry Smith 7583501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 7590452661fSBarry Smith PetscFree(a->rowners); 7603501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A);CHKERRQ(ierr); 7613501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 7623501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 7630452661fSBarry Smith PetscFree(a); 764f830108cSBarry Smith 765f830108cSBarry Smith /* 766f830108cSBarry Smith This is horrible, horrible code. We need to keep the 767f830108cSBarry Smith A pointers for the bops and ops but copy everything 768f830108cSBarry Smith else from C. 769f830108cSBarry Smith */ 770f830108cSBarry Smith Abops = A->bops; 771f830108cSBarry Smith Aops = A->ops; 772*549d3d68SSatish Balay ierr = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr); 773f830108cSBarry Smith A->bops = Abops; 774f830108cSBarry Smith A->ops = Aops; 775f830108cSBarry Smith 7760452661fSBarry Smith PetscHeaderDestroy(B); 7773501a2bdSLois Curfman McInnes } 7783a40ed3dSBarry Smith PetscFunctionReturn(0); 779096963f5SLois Curfman McInnes } 780096963f5SLois Curfman McInnes 781eadb2fb4SBarry Smith #include "pinclude/blaslapack.h" 7825615d1e5SSatish Balay #undef __FUNC__ 7835615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense" 7848f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA) 78544cd7ae7SLois Curfman McInnes { 78644cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 78744cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 78844cd7ae7SLois Curfman McInnes int one = 1, nz; 78944cd7ae7SLois Curfman McInnes 7903a40ed3dSBarry Smith PetscFunctionBegin; 79144cd7ae7SLois Curfman McInnes nz = a->m*a->n; 79244cd7ae7SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 79344cd7ae7SLois Curfman McInnes PLogFlops(nz); 7943a40ed3dSBarry Smith PetscFunctionReturn(0); 79544cd7ae7SLois Curfman McInnes } 79644cd7ae7SLois Curfman McInnes 7975609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 7987b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **); 7998965ea79SLois Curfman McInnes 8008965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 80109dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 80209dc0095SBarry Smith MatGetRow_MPIDense, 80309dc0095SBarry Smith MatRestoreRow_MPIDense, 80409dc0095SBarry Smith MatMult_MPIDense, 80509dc0095SBarry Smith MatMultAdd_MPIDense, 80609dc0095SBarry Smith MatMultTrans_MPIDense, 80709dc0095SBarry Smith MatMultTransAdd_MPIDense, 8088965ea79SLois Curfman McInnes 0, 80909dc0095SBarry Smith 0, 81009dc0095SBarry Smith 0, 81109dc0095SBarry Smith 0, 81209dc0095SBarry Smith 0, 81309dc0095SBarry Smith 0, 81409dc0095SBarry Smith 0, 81509dc0095SBarry Smith MatTranspose_MPIDense, 81609dc0095SBarry Smith MatGetInfo_MPIDense,0, 81709dc0095SBarry Smith MatGetDiagonal_MPIDense, 81809dc0095SBarry Smith 0, 81909dc0095SBarry Smith MatNorm_MPIDense, 82009dc0095SBarry Smith MatAssemblyBegin_MPIDense, 82109dc0095SBarry Smith MatAssemblyEnd_MPIDense, 82209dc0095SBarry Smith 0, 82309dc0095SBarry Smith MatSetOption_MPIDense, 82409dc0095SBarry Smith MatZeroEntries_MPIDense, 82509dc0095SBarry Smith MatZeroRows_MPIDense, 82609dc0095SBarry Smith 0, 82709dc0095SBarry Smith 0, 82809dc0095SBarry Smith 0, 82909dc0095SBarry Smith 0, 83009dc0095SBarry Smith MatGetSize_MPIDense, 83109dc0095SBarry Smith MatGetLocalSize_MPIDense, 83239ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 83309dc0095SBarry Smith 0, 83409dc0095SBarry Smith 0, 83509dc0095SBarry Smith MatGetArray_MPIDense, 83609dc0095SBarry Smith MatRestoreArray_MPIDense, 8375609ef8eSBarry Smith MatDuplicate_MPIDense, 83809dc0095SBarry Smith 0, 83909dc0095SBarry Smith 0, 84009dc0095SBarry Smith 0, 84109dc0095SBarry Smith 0, 84209dc0095SBarry Smith 0, 8432ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 84409dc0095SBarry Smith 0, 84509dc0095SBarry Smith MatGetValues_MPIDense, 84609dc0095SBarry Smith 0, 84709dc0095SBarry Smith 0, 84809dc0095SBarry Smith MatScale_MPIDense, 84909dc0095SBarry Smith 0, 85009dc0095SBarry Smith 0, 85109dc0095SBarry Smith 0, 85209dc0095SBarry Smith MatGetBlockSize_MPIDense, 85309dc0095SBarry Smith 0, 85409dc0095SBarry Smith 0, 85509dc0095SBarry Smith 0, 85609dc0095SBarry Smith 0, 85709dc0095SBarry Smith 0, 85809dc0095SBarry Smith 0, 85909dc0095SBarry Smith 0, 86009dc0095SBarry Smith 0, 86109dc0095SBarry Smith 0, 86209dc0095SBarry Smith 0, 86309dc0095SBarry Smith 0, 86409dc0095SBarry Smith 0, 86509dc0095SBarry Smith MatGetMaps_Petsc}; 8668965ea79SLois Curfman McInnes 8675615d1e5SSatish Balay #undef __FUNC__ 8685615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense" 8698965ea79SLois Curfman McInnes /*@C 87039ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 8718965ea79SLois Curfman McInnes 872db81eaa0SLois Curfman McInnes Collective on MPI_Comm 873db81eaa0SLois Curfman McInnes 8748965ea79SLois Curfman McInnes Input Parameters: 875db81eaa0SLois Curfman McInnes + comm - MPI communicator 8768965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 877db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 8788965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 879db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 880db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 881dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 8828965ea79SLois Curfman McInnes 8838965ea79SLois Curfman McInnes Output Parameter: 884477f1c0bSLois Curfman McInnes . A - the matrix 8858965ea79SLois Curfman McInnes 886b259b22eSLois Curfman McInnes Notes: 88739ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 88839ddd567SLois Curfman McInnes storage by columns. 8898965ea79SLois Curfman McInnes 89018f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 89118f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 892b4fd4287SBarry Smith set data=PETSC_NULL. 89318f449edSLois Curfman McInnes 8948965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 8958965ea79SLois Curfman McInnes (possibly both). 8968965ea79SLois Curfman McInnes 8973501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 8983501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 8993501a2bdSLois Curfman McInnes 900027ccd11SLois Curfman McInnes Level: intermediate 901027ccd11SLois Curfman McInnes 90239ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 9038965ea79SLois Curfman McInnes 90439ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 9058965ea79SLois Curfman McInnes @*/ 906477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 9078965ea79SLois Curfman McInnes { 9088965ea79SLois Curfman McInnes Mat mat; 90939ddd567SLois Curfman McInnes Mat_MPIDense *a; 91025cdf11fSBarry Smith int ierr, i,flg; 9118965ea79SLois Curfman McInnes 9123a40ed3dSBarry Smith PetscFunctionBegin; 913ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 914ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 91518f449edSLois Curfman McInnes 916477f1c0bSLois Curfman McInnes *A = 0; 9173f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView); 9188965ea79SLois Curfman McInnes PLogObjectCreate(mat); 9190452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 920*549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 921e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 922e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 9238965ea79SLois Curfman McInnes mat->factor = 0; 92490f02eecSBarry Smith mat->mapping = 0; 9258965ea79SLois Curfman McInnes 926622d7880SLois Curfman McInnes a->factor = 0; 927e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 928d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr); 929d132466eSBarry Smith ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr); 9308965ea79SLois Curfman McInnes 93196f6c058SBarry Smith ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr); 93239ddd567SLois Curfman McInnes 933c7fcc2eaSBarry Smith /* 934c7fcc2eaSBarry Smith The computation of n is wrong below, n should represent the number of local 935c7fcc2eaSBarry Smith rows in the right (column vector) 936c7fcc2eaSBarry Smith */ 937c7fcc2eaSBarry Smith 93839ddd567SLois Curfman McInnes /* each row stores all columns */ 93939ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 940d7e8b826SBarry Smith if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 941a8c6a408SBarry Smith /* if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */ 942aca0ad90SLois Curfman McInnes a->N = mat->N = N; 943aca0ad90SLois Curfman McInnes a->M = mat->M = M; 944aca0ad90SLois Curfman McInnes a->m = mat->m = m; 945aca0ad90SLois Curfman McInnes a->n = mat->n = n; 9468965ea79SLois Curfman McInnes 947c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 948c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 949488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr); 95096f6c058SBarry Smith ierr = MapCreateMPI(comm,PETSC_DECIDE,N,&mat->cmap);CHKERRQ(ierr); 951c7fcc2eaSBarry Smith 9528965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 953d7e8b826SBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners); 954d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 955f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 956ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 9578965ea79SLois Curfman McInnes a->rowners[0] = 0; 9588965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 9598965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 9608965ea79SLois Curfman McInnes } 9618965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 9628965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 963ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 964d7e8b826SBarry Smith a->cowners[0] = 0; 965d7e8b826SBarry Smith for ( i=2; i<=a->size; i++ ) { 966d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 967d7e8b826SBarry Smith } 9688965ea79SLois Curfman McInnes 969029af93fSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr); 9708965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 9718965ea79SLois Curfman McInnes 9728965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 9733782ba37SSatish Balay a->donotstash = 0; 9748798bf22SSatish Balay ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr); 9758965ea79SLois Curfman McInnes 9768965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 9778965ea79SLois Curfman McInnes a->lvec = 0; 9788965ea79SLois Curfman McInnes a->Mvctx = 0; 97939b7565bSBarry Smith a->roworiented = 1; 9808965ea79SLois Curfman McInnes 981477f1c0bSLois Curfman McInnes *A = mat; 98225cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 98325cdf11fSBarry Smith if (flg) { 9848c469469SLois Curfman McInnes ierr = MatPrintHelp(mat);CHKERRQ(ierr); 9858c469469SLois Curfman McInnes } 9863a40ed3dSBarry Smith PetscFunctionReturn(0); 9878965ea79SLois Curfman McInnes } 9888965ea79SLois Curfman McInnes 9895615d1e5SSatish Balay #undef __FUNC__ 9905609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense" 9915609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9928965ea79SLois Curfman McInnes { 9938965ea79SLois Curfman McInnes Mat mat; 9943501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 99539ddd567SLois Curfman McInnes int ierr; 9962ba99913SLois Curfman McInnes FactorCtx *factor; 9978965ea79SLois Curfman McInnes 9983a40ed3dSBarry Smith PetscFunctionBegin; 9998965ea79SLois Curfman McInnes *newmat = 0; 10003f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView); 10018965ea79SLois Curfman McInnes PLogObjectCreate(mat); 10020452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 1003*549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1004e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1005e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 10063501a2bdSLois Curfman McInnes mat->factor = A->factor; 1007c456f294SBarry Smith mat->assembled = PETSC_TRUE; 10088965ea79SLois Curfman McInnes 100944cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 101044cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 101144cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 101244cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 10132ba99913SLois Curfman McInnes if (oldmat->factor) { 10142ba99913SLois Curfman McInnes a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx));CHKPTRQ(factor); 10152ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 10162ba99913SLois Curfman McInnes } else a->factor = 0; 10178965ea79SLois Curfman McInnes 10188965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 10198965ea79SLois Curfman McInnes a->rend = oldmat->rend; 10208965ea79SLois Curfman McInnes a->size = oldmat->size; 10218965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1022e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 10233782ba37SSatish Balay a->donotstash = oldmat->donotstash; 10240452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners); 1025f09e8eb9SSatish Balay PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1026*549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 10278798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 10288965ea79SLois Curfman McInnes 10298965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr); 10308965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 103155659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr); 10328965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 10335609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 10348965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10358965ea79SLois Curfman McInnes *newmat = mat; 10363a40ed3dSBarry Smith PetscFunctionReturn(0); 10378965ea79SLois Curfman McInnes } 10388965ea79SLois Curfman McInnes 103977c4ece6SBarry Smith #include "sys.h" 10408965ea79SLois Curfman McInnes 10415615d1e5SSatish Balay #undef __FUNC__ 10425615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile" 104390ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 104490ace30eSBarry Smith { 104540011551SBarry Smith int *rowners, i,size,rank,m,ierr,nz,j; 104690ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 104790ace30eSBarry Smith MPI_Status status; 104890ace30eSBarry Smith 10493a40ed3dSBarry Smith PetscFunctionBegin; 1050d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1051d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 105290ace30eSBarry Smith 105390ace30eSBarry Smith /* determine ownership of all rows */ 105490ace30eSBarry Smith m = M/size + ((M % size) > rank); 105590ace30eSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1056ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 105790ace30eSBarry Smith rowners[0] = 0; 105890ace30eSBarry Smith for ( i=2; i<=size; i++ ) { 105990ace30eSBarry Smith rowners[i] += rowners[i-1]; 106090ace30eSBarry Smith } 106190ace30eSBarry Smith 106290ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 106390ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 106490ace30eSBarry Smith 106590ace30eSBarry Smith if (!rank) { 106690ace30eSBarry Smith vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 106790ace30eSBarry Smith 106890ace30eSBarry Smith /* read in my part of the matrix numerical values */ 10690752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 107090ace30eSBarry Smith 107190ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 107290ace30eSBarry Smith vals_ptr = vals; 107390ace30eSBarry Smith for ( i=0; i<m; i++ ) { 107490ace30eSBarry Smith for ( j=0; j<N; j++ ) { 107590ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 107690ace30eSBarry Smith } 107790ace30eSBarry Smith } 107890ace30eSBarry Smith 107990ace30eSBarry Smith /* read in other processors and ship out */ 108090ace30eSBarry Smith for ( i=1; i<size; i++ ) { 108190ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 10820752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1083ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 108490ace30eSBarry Smith } 10853a40ed3dSBarry Smith } else { 108690ace30eSBarry Smith /* receive numeric values */ 108790ace30eSBarry Smith vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 108890ace30eSBarry Smith 108990ace30eSBarry Smith /* receive message of values*/ 1090ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 109190ace30eSBarry Smith 109290ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 109390ace30eSBarry Smith vals_ptr = vals; 109490ace30eSBarry Smith for ( i=0; i<m; i++ ) { 109590ace30eSBarry Smith for ( j=0; j<N; j++ ) { 109690ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 109790ace30eSBarry Smith } 109890ace30eSBarry Smith } 109990ace30eSBarry Smith } 110090ace30eSBarry Smith PetscFree(rowners); 110190ace30eSBarry Smith PetscFree(vals); 11026d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11036d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11043a40ed3dSBarry Smith PetscFunctionReturn(0); 110590ace30eSBarry Smith } 110690ace30eSBarry Smith 110790ace30eSBarry Smith 11085615d1e5SSatish Balay #undef __FUNC__ 11095615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense" 111019bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 11118965ea79SLois Curfman McInnes { 11128965ea79SLois Curfman McInnes Mat A; 11138965ea79SLois Curfman McInnes Scalar *vals,*svals; 111419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 11158965ea79SLois Curfman McInnes MPI_Status status; 11168965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 11178965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 111819bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 11193a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 11208965ea79SLois Curfman McInnes 11213a40ed3dSBarry Smith PetscFunctionBegin; 1122d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1123d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 11248965ea79SLois Curfman McInnes if (!rank) { 112590ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 11260752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1127a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 11288965ea79SLois Curfman McInnes } 11298965ea79SLois Curfman McInnes 1130ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 113190ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 113290ace30eSBarry Smith 113390ace30eSBarry Smith /* 113490ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 113590ace30eSBarry Smith */ 113690ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 11373a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 11383a40ed3dSBarry Smith PetscFunctionReturn(0); 113990ace30eSBarry Smith } 114090ace30eSBarry Smith 11418965ea79SLois Curfman McInnes /* determine ownership of all rows */ 11428965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 11430452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1144ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 11458965ea79SLois Curfman McInnes rowners[0] = 0; 11468965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 11478965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 11488965ea79SLois Curfman McInnes } 11498965ea79SLois Curfman McInnes rstart = rowners[rank]; 11508965ea79SLois Curfman McInnes rend = rowners[rank+1]; 11518965ea79SLois Curfman McInnes 11528965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 11530452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) );CHKPTRQ(ourlens); 11548965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 11558965ea79SLois Curfman McInnes if (!rank) { 11560452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) );CHKPTRQ(rowlengths); 11570752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 11580452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(sndcounts); 11598965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1160ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 11610452661fSBarry Smith PetscFree(sndcounts); 1162ca161407SBarry Smith } else { 1163ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 11648965ea79SLois Curfman McInnes } 11658965ea79SLois Curfman McInnes 11668965ea79SLois Curfman McInnes if (!rank) { 11678965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 11680452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(procsnz); 1169*549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 11708965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 11718965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 11728965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 11738965ea79SLois Curfman McInnes } 11748965ea79SLois Curfman McInnes } 11750452661fSBarry Smith PetscFree(rowlengths); 11768965ea79SLois Curfman McInnes 11778965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 11788965ea79SLois Curfman McInnes maxnz = 0; 11798965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 11800452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 11818965ea79SLois Curfman McInnes } 11820452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) );CHKPTRQ(cols); 11838965ea79SLois Curfman McInnes 11848965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 11858965ea79SLois Curfman McInnes nz = procsnz[0]; 11860452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 11870752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 11888965ea79SLois Curfman McInnes 11898965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 11908965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 11918965ea79SLois Curfman McInnes nz = procsnz[i]; 11920752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1193ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 11948965ea79SLois Curfman McInnes } 11950452661fSBarry Smith PetscFree(cols); 11963a40ed3dSBarry Smith } else { 11978965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 11988965ea79SLois Curfman McInnes nz = 0; 11998965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12008965ea79SLois Curfman McInnes nz += ourlens[i]; 12018965ea79SLois Curfman McInnes } 12020452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 12038965ea79SLois Curfman McInnes 12048965ea79SLois Curfman McInnes /* receive message of column indices*/ 1205ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1206ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1207a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12088965ea79SLois Curfman McInnes } 12098965ea79SLois Curfman McInnes 12108965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1211*549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 12128965ea79SLois Curfman McInnes jj = 0; 12138965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12148965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 12158965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 12168965ea79SLois Curfman McInnes jj++; 12178965ea79SLois Curfman McInnes } 12188965ea79SLois Curfman McInnes } 12198965ea79SLois Curfman McInnes 12208965ea79SLois Curfman McInnes /* create our matrix */ 12218965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12228965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 12238965ea79SLois Curfman McInnes } 1224b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 12258965ea79SLois Curfman McInnes A = *newmat; 12268965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12278965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 12288965ea79SLois Curfman McInnes } 12298965ea79SLois Curfman McInnes 12308965ea79SLois Curfman McInnes if (!rank) { 12310452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) );CHKPTRQ(vals); 12328965ea79SLois Curfman McInnes 12338965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 12348965ea79SLois Curfman McInnes nz = procsnz[0]; 12350752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 12368965ea79SLois Curfman McInnes 12378965ea79SLois Curfman McInnes /* insert into matrix */ 12388965ea79SLois Curfman McInnes jj = rstart; 12398965ea79SLois Curfman McInnes smycols = mycols; 12408965ea79SLois Curfman McInnes svals = vals; 12418965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12428965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12438965ea79SLois Curfman McInnes smycols += ourlens[i]; 12448965ea79SLois Curfman McInnes svals += ourlens[i]; 12458965ea79SLois Curfman McInnes jj++; 12468965ea79SLois Curfman McInnes } 12478965ea79SLois Curfman McInnes 12488965ea79SLois Curfman McInnes /* read in other processors and ship out */ 12498965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 12508965ea79SLois Curfman McInnes nz = procsnz[i]; 12510752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1252ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 12538965ea79SLois Curfman McInnes } 12540452661fSBarry Smith PetscFree(procsnz); 12553a40ed3dSBarry Smith } else { 12568965ea79SLois Curfman McInnes /* receive numeric values */ 12570452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) );CHKPTRQ(vals); 12588965ea79SLois Curfman McInnes 12598965ea79SLois Curfman McInnes /* receive message of values*/ 1260ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1261ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1262a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12638965ea79SLois Curfman McInnes 12648965ea79SLois Curfman McInnes /* insert into matrix */ 12658965ea79SLois Curfman McInnes jj = rstart; 12668965ea79SLois Curfman McInnes smycols = mycols; 12678965ea79SLois Curfman McInnes svals = vals; 12688965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12698965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12708965ea79SLois Curfman McInnes smycols += ourlens[i]; 12718965ea79SLois Curfman McInnes svals += ourlens[i]; 12728965ea79SLois Curfman McInnes jj++; 12738965ea79SLois Curfman McInnes } 12748965ea79SLois Curfman McInnes } 12750452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 12768965ea79SLois Curfman McInnes 12776d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12786d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12793a40ed3dSBarry Smith PetscFunctionReturn(0); 12808965ea79SLois Curfman McInnes } 128190ace30eSBarry Smith 128290ace30eSBarry Smith 128390ace30eSBarry Smith 128490ace30eSBarry Smith 128590ace30eSBarry Smith 1286